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 18413f72f0SBarry Smith 19*9305a4c7SMatthew G. Knepley PetscErrorCode ISGetPointRange(IS pointIS, PetscInt *pStart, PetscInt *pEnd, const PetscInt **points) 20*9305a4c7SMatthew G. Knepley { 21*9305a4c7SMatthew G. Knepley PetscInt numCells, step = 1; 22*9305a4c7SMatthew G. Knepley PetscBool isStride; 23*9305a4c7SMatthew G. Knepley PetscErrorCode ierr; 24*9305a4c7SMatthew G. Knepley 25*9305a4c7SMatthew G. Knepley PetscFunctionBeginHot; 26*9305a4c7SMatthew G. Knepley *pStart = 0; 27*9305a4c7SMatthew G. Knepley *points = NULL; 28*9305a4c7SMatthew G. Knepley ierr = ISGetLocalSize(pointIS, &numCells);CHKERRQ(ierr); 29*9305a4c7SMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject) pointIS, ISSTRIDE, &isStride);CHKERRQ(ierr); 30*9305a4c7SMatthew G. Knepley if (isStride) {ierr = ISStrideGetInfo(pointIS, pStart, &step);CHKERRQ(ierr);} 31*9305a4c7SMatthew G. Knepley *pEnd = *pStart + numCells; 32*9305a4c7SMatthew G. Knepley if (!isStride || step != 1) {ierr = ISGetIndices(pointIS, points);CHKERRQ(ierr);} 33*9305a4c7SMatthew G. Knepley PetscFunctionReturn(0); 34*9305a4c7SMatthew G. Knepley } 35*9305a4c7SMatthew G. Knepley 36*9305a4c7SMatthew G. Knepley PetscErrorCode ISRestorePointRange(IS pointIS, PetscInt *pStart, PetscInt *pEnd, const PetscInt **points) 37*9305a4c7SMatthew G. Knepley { 38*9305a4c7SMatthew G. Knepley PetscInt step = 1; 39*9305a4c7SMatthew G. Knepley PetscBool isStride; 40*9305a4c7SMatthew G. Knepley PetscErrorCode ierr; 41*9305a4c7SMatthew G. Knepley 42*9305a4c7SMatthew G. Knepley PetscFunctionBeginHot; 43*9305a4c7SMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject) pointIS, ISSTRIDE, &isStride);CHKERRQ(ierr); 44*9305a4c7SMatthew G. Knepley if (isStride) {ierr = ISStrideGetInfo(pointIS, pStart, &step);CHKERRQ(ierr);} 45*9305a4c7SMatthew G. Knepley if (!isStride || step != 1) {ierr = ISGetIndices(pointIS, points);CHKERRQ(ierr);} 46*9305a4c7SMatthew G. Knepley PetscFunctionReturn(0); 47*9305a4c7SMatthew G. Knepley } 48*9305a4c7SMatthew G. Knepley 49*9305a4c7SMatthew G. Knepley PetscErrorCode ISGetPointSubrange(IS subpointIS, PetscInt pStart, PetscInt pEnd, const PetscInt *points) 50*9305a4c7SMatthew G. Knepley { 51*9305a4c7SMatthew G. Knepley PetscErrorCode ierr; 52*9305a4c7SMatthew G. Knepley 53*9305a4c7SMatthew G. Knepley PetscFunctionBeginHot; 54*9305a4c7SMatthew G. Knepley if (points) { 55*9305a4c7SMatthew G. Knepley ierr = ISSetType(subpointIS, ISGENERAL);CHKERRQ(ierr); 56*9305a4c7SMatthew G. Knepley ierr = ISGeneralSetIndices(subpointIS, pEnd-pStart, &points[pStart], PETSC_USE_POINTER);CHKERRQ(ierr); 57*9305a4c7SMatthew G. Knepley } else { 58*9305a4c7SMatthew G. Knepley ierr = ISSetType(subpointIS, ISSTRIDE);CHKERRQ(ierr); 59*9305a4c7SMatthew G. Knepley ierr = ISStrideSetStride(subpointIS, pEnd-pStart, pStart, 1);CHKERRQ(ierr); 60*9305a4c7SMatthew G. Knepley } 61*9305a4c7SMatthew G. Knepley PetscFunctionReturn(0); 62*9305a4c7SMatthew G. Knepley } 63*9305a4c7SMatthew G. Knepley 64413f72f0SBarry Smith /* -----------------------------------------------------------------------------------------*/ 65413f72f0SBarry Smith 66413f72f0SBarry Smith /* 67413f72f0SBarry Smith Creates the global mapping information in the ISLocalToGlobalMapping structure 68413f72f0SBarry Smith 69413f72f0SBarry Smith If the user has not selected how to handle the global to local mapping then use HASH for "large" problems 70413f72f0SBarry Smith */ 71413f72f0SBarry Smith static PetscErrorCode ISGlobalToLocalMappingSetUp(ISLocalToGlobalMapping mapping) 72413f72f0SBarry Smith { 73413f72f0SBarry Smith PetscInt i,*idx = mapping->indices,n = mapping->n,end,start; 74413f72f0SBarry Smith PetscErrorCode ierr; 75413f72f0SBarry Smith 76413f72f0SBarry Smith PetscFunctionBegin; 77413f72f0SBarry Smith if (mapping->data) PetscFunctionReturn(0); 78413f72f0SBarry Smith end = 0; 79413f72f0SBarry Smith start = PETSC_MAX_INT; 80413f72f0SBarry Smith 81413f72f0SBarry Smith for (i=0; i<n; i++) { 82413f72f0SBarry Smith if (idx[i] < 0) continue; 83413f72f0SBarry Smith if (idx[i] < start) start = idx[i]; 84413f72f0SBarry Smith if (idx[i] > end) end = idx[i]; 85413f72f0SBarry Smith } 86413f72f0SBarry Smith if (start > end) {start = 0; end = -1;} 87413f72f0SBarry Smith mapping->globalstart = start; 88413f72f0SBarry Smith mapping->globalend = end; 89413f72f0SBarry Smith if (!((PetscObject)mapping)->type_name) { 90413f72f0SBarry Smith if ((end - start) > PetscMax(4*n,1000000)) { 917f79407eSBarry Smith ierr = ISLocalToGlobalMappingSetType(mapping,ISLOCALTOGLOBALMAPPINGHASH);CHKERRQ(ierr); 92413f72f0SBarry Smith } else { 937f79407eSBarry Smith ierr = ISLocalToGlobalMappingSetType(mapping,ISLOCALTOGLOBALMAPPINGBASIC);CHKERRQ(ierr); 94413f72f0SBarry Smith } 95413f72f0SBarry Smith } 96413f72f0SBarry Smith ierr = (*mapping->ops->globaltolocalmappingsetup)(mapping);CHKERRQ(ierr); 97413f72f0SBarry Smith PetscFunctionReturn(0); 98413f72f0SBarry Smith } 99413f72f0SBarry Smith 100413f72f0SBarry Smith static PetscErrorCode ISGlobalToLocalMappingSetUp_Basic(ISLocalToGlobalMapping mapping) 101413f72f0SBarry Smith { 102413f72f0SBarry Smith PetscErrorCode ierr; 103413f72f0SBarry Smith PetscInt i,*idx = mapping->indices,n = mapping->n,end,start,*globals; 104413f72f0SBarry Smith ISLocalToGlobalMapping_Basic *map; 105413f72f0SBarry Smith 106413f72f0SBarry Smith PetscFunctionBegin; 107413f72f0SBarry Smith start = mapping->globalstart; 108413f72f0SBarry Smith end = mapping->globalend; 109413f72f0SBarry Smith ierr = PetscNew(&map);CHKERRQ(ierr); 110413f72f0SBarry Smith ierr = PetscMalloc1(end-start+2,&globals);CHKERRQ(ierr); 111413f72f0SBarry Smith map->globals = globals; 112413f72f0SBarry Smith for (i=0; i<end-start+1; i++) globals[i] = -1; 113413f72f0SBarry Smith for (i=0; i<n; i++) { 114413f72f0SBarry Smith if (idx[i] < 0) continue; 115413f72f0SBarry Smith globals[idx[i] - start] = i; 116413f72f0SBarry Smith } 117413f72f0SBarry Smith mapping->data = (void*)map; 118413f72f0SBarry Smith ierr = PetscLogObjectMemory((PetscObject)mapping,(end-start+1)*sizeof(PetscInt));CHKERRQ(ierr); 119413f72f0SBarry Smith PetscFunctionReturn(0); 120413f72f0SBarry Smith } 121413f72f0SBarry Smith 122413f72f0SBarry Smith static PetscErrorCode ISGlobalToLocalMappingSetUp_Hash(ISLocalToGlobalMapping mapping) 123413f72f0SBarry Smith { 124413f72f0SBarry Smith PetscErrorCode ierr; 125413f72f0SBarry Smith PetscInt i,*idx = mapping->indices,n = mapping->n; 126413f72f0SBarry Smith ISLocalToGlobalMapping_Hash *map; 127413f72f0SBarry Smith 128413f72f0SBarry Smith PetscFunctionBegin; 129413f72f0SBarry Smith ierr = PetscNew(&map);CHKERRQ(ierr); 130e8f14785SLisandro Dalcin ierr = PetscHMapICreate(&map->globalht);CHKERRQ(ierr); 131413f72f0SBarry Smith for (i=0; i<n; i++ ) { 132413f72f0SBarry Smith if (idx[i] < 0) continue; 133e8f14785SLisandro Dalcin ierr = PetscHMapISet(map->globalht,idx[i],i);CHKERRQ(ierr); 134413f72f0SBarry Smith } 135413f72f0SBarry Smith mapping->data = (void*)map; 136413f72f0SBarry Smith ierr = PetscLogObjectMemory((PetscObject)mapping,2*n*sizeof(PetscInt));CHKERRQ(ierr); 137413f72f0SBarry Smith PetscFunctionReturn(0); 138413f72f0SBarry Smith } 139413f72f0SBarry Smith 140413f72f0SBarry Smith static PetscErrorCode ISLocalToGlobalMappingDestroy_Basic(ISLocalToGlobalMapping mapping) 141413f72f0SBarry Smith { 142413f72f0SBarry Smith PetscErrorCode ierr; 143413f72f0SBarry Smith ISLocalToGlobalMapping_Basic *map = (ISLocalToGlobalMapping_Basic *)mapping->data; 144413f72f0SBarry Smith 145413f72f0SBarry Smith PetscFunctionBegin; 146413f72f0SBarry Smith if (!map) PetscFunctionReturn(0); 147413f72f0SBarry Smith ierr = PetscFree(map->globals);CHKERRQ(ierr); 148413f72f0SBarry Smith ierr = PetscFree(mapping->data);CHKERRQ(ierr); 149413f72f0SBarry Smith PetscFunctionReturn(0); 150413f72f0SBarry Smith } 151413f72f0SBarry Smith 152413f72f0SBarry Smith static PetscErrorCode ISLocalToGlobalMappingDestroy_Hash(ISLocalToGlobalMapping mapping) 153413f72f0SBarry Smith { 154413f72f0SBarry Smith PetscErrorCode ierr; 155413f72f0SBarry Smith ISLocalToGlobalMapping_Hash *map = (ISLocalToGlobalMapping_Hash*)mapping->data; 156413f72f0SBarry Smith 157413f72f0SBarry Smith PetscFunctionBegin; 158413f72f0SBarry Smith if (!map) PetscFunctionReturn(0); 159e8f14785SLisandro Dalcin ierr = PetscHMapIDestroy(&map->globalht);CHKERRQ(ierr); 160413f72f0SBarry Smith ierr = PetscFree(mapping->data);CHKERRQ(ierr); 161413f72f0SBarry Smith PetscFunctionReturn(0); 162413f72f0SBarry Smith } 163413f72f0SBarry Smith 164413f72f0SBarry Smith #define GTOLTYPE _Basic 165413f72f0SBarry Smith #define GTOLNAME _Basic 166541bf97eSAdrian Croucher #define GTOLBS mapping->bs 167413f72f0SBarry Smith #define GTOL(g, local) do { \ 168413f72f0SBarry Smith local = map->globals[g/bs - start]; \ 169413f72f0SBarry Smith local = bs*local + (g % bs); \ 170413f72f0SBarry Smith } while (0) 171541bf97eSAdrian Croucher 172413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h> 173413f72f0SBarry Smith 174413f72f0SBarry Smith #define GTOLTYPE _Basic 175413f72f0SBarry Smith #define GTOLNAME Block_Basic 176541bf97eSAdrian Croucher #define GTOLBS 1 177413f72f0SBarry Smith #define GTOL(g, local) do { \ 178413f72f0SBarry Smith local = map->globals[g - start]; \ 179413f72f0SBarry Smith } while (0) 180413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h> 181413f72f0SBarry Smith 182413f72f0SBarry Smith #define GTOLTYPE _Hash 183413f72f0SBarry Smith #define GTOLNAME _Hash 184541bf97eSAdrian Croucher #define GTOLBS mapping->bs 185413f72f0SBarry Smith #define GTOL(g, local) do { \ 186e8f14785SLisandro Dalcin (void)PetscHMapIGet(map->globalht,g/bs,&local); \ 187413f72f0SBarry Smith local = bs*local + (g % bs); \ 188413f72f0SBarry Smith } while (0) 189413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h> 190413f72f0SBarry Smith 191413f72f0SBarry Smith #define GTOLTYPE _Hash 192413f72f0SBarry Smith #define GTOLNAME Block_Hash 193541bf97eSAdrian Croucher #define GTOLBS 1 194413f72f0SBarry Smith #define GTOL(g, local) do { \ 195e8f14785SLisandro Dalcin (void)PetscHMapIGet(map->globalht,g,&local); \ 196413f72f0SBarry Smith } while (0) 197413f72f0SBarry Smith #include <../src/vec/is/utils/isltog.h> 198413f72f0SBarry Smith 1996658fb44Sstefano_zampini /*@ 2006658fb44Sstefano_zampini ISLocalToGlobalMappingDuplicate - Duplicates the local to global mapping object 2016658fb44Sstefano_zampini 2026658fb44Sstefano_zampini Not Collective 2036658fb44Sstefano_zampini 2046658fb44Sstefano_zampini Input Parameter: 2056658fb44Sstefano_zampini . ltog - local to global mapping 2066658fb44Sstefano_zampini 2076658fb44Sstefano_zampini Output Parameter: 2086658fb44Sstefano_zampini . nltog - the duplicated local to global mapping 2096658fb44Sstefano_zampini 2106658fb44Sstefano_zampini Level: advanced 2116658fb44Sstefano_zampini 2126658fb44Sstefano_zampini Concepts: mapping^local to global 2136658fb44Sstefano_zampini 2146658fb44Sstefano_zampini .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreate() 2156658fb44Sstefano_zampini @*/ 2166658fb44Sstefano_zampini PetscErrorCode ISLocalToGlobalMappingDuplicate(ISLocalToGlobalMapping ltog,ISLocalToGlobalMapping* nltog) 2176658fb44Sstefano_zampini { 2186658fb44Sstefano_zampini PetscErrorCode ierr; 2196658fb44Sstefano_zampini 2206658fb44Sstefano_zampini PetscFunctionBegin; 2216658fb44Sstefano_zampini PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1); 2226658fb44Sstefano_zampini ierr = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)ltog),ltog->bs,ltog->n,ltog->indices,PETSC_COPY_VALUES,nltog);CHKERRQ(ierr); 2236658fb44Sstefano_zampini PetscFunctionReturn(0); 2246658fb44Sstefano_zampini } 2256658fb44Sstefano_zampini 226565245c5SBarry Smith /*@ 227107e9a97SBarry Smith ISLocalToGlobalMappingGetSize - Gets the local size of a local to global mapping 2283b9aefa3SBarry Smith 2293b9aefa3SBarry Smith Not Collective 2303b9aefa3SBarry Smith 2313b9aefa3SBarry Smith Input Parameter: 2323b9aefa3SBarry Smith . ltog - local to global mapping 2333b9aefa3SBarry Smith 2343b9aefa3SBarry Smith Output Parameter: 235107e9a97SBarry Smith . n - the number of entries in the local mapping, ISLocalToGlobalMappingGetIndices() returns an array of this length 2363b9aefa3SBarry Smith 2373b9aefa3SBarry Smith Level: advanced 2383b9aefa3SBarry Smith 239273d9f13SBarry Smith Concepts: mapping^local to global 2403b9aefa3SBarry Smith 2413b9aefa3SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreate() 2423b9aefa3SBarry Smith @*/ 2437087cfbeSBarry Smith PetscErrorCode ISLocalToGlobalMappingGetSize(ISLocalToGlobalMapping mapping,PetscInt *n) 2443b9aefa3SBarry Smith { 2453b9aefa3SBarry Smith PetscFunctionBegin; 2460700a824SBarry Smith PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1); 2474482741eSBarry Smith PetscValidIntPointer(n,2); 248107e9a97SBarry Smith *n = mapping->bs*mapping->n; 2493b9aefa3SBarry Smith PetscFunctionReturn(0); 2503b9aefa3SBarry Smith } 2513b9aefa3SBarry Smith 2525a5d4f66SBarry Smith /*@C 2535a5d4f66SBarry Smith ISLocalToGlobalMappingView - View a local to global mapping 2545a5d4f66SBarry Smith 255b9cd556bSLois Curfman McInnes Not Collective 256b9cd556bSLois Curfman McInnes 2575a5d4f66SBarry Smith Input Parameters: 2583b9aefa3SBarry Smith + ltog - local to global mapping 2593b9aefa3SBarry Smith - viewer - viewer 2605a5d4f66SBarry Smith 261a997ad1aSLois Curfman McInnes Level: advanced 262a997ad1aSLois Curfman McInnes 263273d9f13SBarry Smith Concepts: mapping^local to global 2645a5d4f66SBarry Smith 2655a5d4f66SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreate() 2665a5d4f66SBarry Smith @*/ 2677087cfbeSBarry Smith PetscErrorCode ISLocalToGlobalMappingView(ISLocalToGlobalMapping mapping,PetscViewer viewer) 2685a5d4f66SBarry Smith { 26932dcc486SBarry Smith PetscInt i; 27032dcc486SBarry Smith PetscMPIInt rank; 271ace3abfcSBarry Smith PetscBool iascii; 2726849ba73SBarry Smith PetscErrorCode ierr; 2735a5d4f66SBarry Smith 2745a5d4f66SBarry Smith PetscFunctionBegin; 2750700a824SBarry Smith PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1); 2763050cee2SBarry Smith if (!viewer) { 277ce94432eSBarry Smith ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)mapping),&viewer);CHKERRQ(ierr); 2783050cee2SBarry Smith } 2790700a824SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 2805a5d4f66SBarry Smith 281ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)mapping),&rank);CHKERRQ(ierr); 282251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 28332077d6dSBarry Smith if (iascii) { 28498c3331eSBarry Smith ierr = PetscObjectPrintClassNamePrefixType((PetscObject)mapping,viewer);CHKERRQ(ierr); 2851575c14dSBarry Smith ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 2865a5d4f66SBarry Smith for (i=0; i<mapping->n; i++) { 2877904a332SBarry Smith ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d] %D %D\n",rank,i,mapping->indices[i]);CHKERRQ(ierr); 2886831982aSBarry Smith } 289b0a32e0cSBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 2901575c14dSBarry Smith ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 2911575c14dSBarry Smith } 2925a5d4f66SBarry Smith PetscFunctionReturn(0); 2935a5d4f66SBarry Smith } 2945a5d4f66SBarry Smith 2951f428162SBarry Smith /*@ 2962bdab257SBarry Smith ISLocalToGlobalMappingCreateIS - Creates a mapping between a local (0 to n) 2972bdab257SBarry Smith ordering and a global parallel ordering. 2982bdab257SBarry Smith 2990f5bd95cSBarry Smith Not collective 300b9cd556bSLois Curfman McInnes 301a997ad1aSLois Curfman McInnes Input Parameter: 3028c03b21aSDmitry Karpeev . is - index set containing the global numbers for each local number 3032bdab257SBarry Smith 304a997ad1aSLois Curfman McInnes Output Parameter: 3052bdab257SBarry Smith . mapping - new mapping data structure 3062bdab257SBarry Smith 30795452b02SPatrick Sanan Notes: 30895452b02SPatrick Sanan the block size of the IS determines the block size of the mapping 309a997ad1aSLois Curfman McInnes Level: advanced 310a997ad1aSLois Curfman McInnes 311273d9f13SBarry Smith Concepts: mapping^local to global 3122bdab257SBarry Smith 3137e99dc12SLawrence Mitchell .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingSetFromOptions() 3142bdab257SBarry Smith @*/ 3157087cfbeSBarry Smith PetscErrorCode ISLocalToGlobalMappingCreateIS(IS is,ISLocalToGlobalMapping *mapping) 3162bdab257SBarry Smith { 3176849ba73SBarry Smith PetscErrorCode ierr; 3183bbf0e92SBarry Smith PetscInt n,bs; 3195d0c19d7SBarry Smith const PetscInt *indices; 3202bdab257SBarry Smith MPI_Comm comm; 3213bbf0e92SBarry Smith PetscBool isblock; 3223a40ed3dSBarry Smith 3233a40ed3dSBarry Smith PetscFunctionBegin; 3240700a824SBarry Smith PetscValidHeaderSpecific(is,IS_CLASSID,1); 3254482741eSBarry Smith PetscValidPointer(mapping,2); 3262bdab257SBarry Smith 3272bdab257SBarry Smith ierr = PetscObjectGetComm((PetscObject)is,&comm);CHKERRQ(ierr); 3283b9aefa3SBarry Smith ierr = ISGetLocalSize(is,&n);CHKERRQ(ierr); 3293bbf0e92SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)is,ISBLOCK,&isblock);CHKERRQ(ierr); 3306006e8d2SBarry Smith if (!isblock) { 331f0413b6fSBarry Smith ierr = ISGetIndices(is,&indices);CHKERRQ(ierr); 332f0413b6fSBarry Smith ierr = ISLocalToGlobalMappingCreate(comm,1,n,indices,PETSC_COPY_VALUES,mapping);CHKERRQ(ierr); 3332bdab257SBarry Smith ierr = ISRestoreIndices(is,&indices);CHKERRQ(ierr); 3346006e8d2SBarry Smith } else { 3356006e8d2SBarry Smith ierr = ISGetBlockSize(is,&bs);CHKERRQ(ierr); 336f0413b6fSBarry Smith ierr = ISBlockGetIndices(is,&indices);CHKERRQ(ierr); 33728bc9809SBarry Smith ierr = ISLocalToGlobalMappingCreate(comm,bs,n/bs,indices,PETSC_COPY_VALUES,mapping);CHKERRQ(ierr); 338f0413b6fSBarry Smith ierr = ISBlockRestoreIndices(is,&indices);CHKERRQ(ierr); 3396006e8d2SBarry Smith } 3403a40ed3dSBarry Smith PetscFunctionReturn(0); 3412bdab257SBarry Smith } 3425a5d4f66SBarry Smith 343a4d96a55SJed Brown /*@C 344a4d96a55SJed Brown ISLocalToGlobalMappingCreateSF - Creates a mapping between a local (0 to n) 345a4d96a55SJed Brown ordering and a global parallel ordering. 346a4d96a55SJed Brown 347a4d96a55SJed Brown Collective 348a4d96a55SJed Brown 349a4d96a55SJed Brown Input Parameter: 350a4d96a55SJed Brown + sf - star forest mapping contiguous local indices to (rank, offset) 351a4d96a55SJed Brown - start - first global index on this process 352a4d96a55SJed Brown 353a4d96a55SJed Brown Output Parameter: 354a4d96a55SJed Brown . mapping - new mapping data structure 355a4d96a55SJed Brown 356a4d96a55SJed Brown Level: advanced 357a4d96a55SJed Brown 358a4d96a55SJed Brown Concepts: mapping^local to global 359a4d96a55SJed Brown 3607e99dc12SLawrence Mitchell .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingSetFromOptions() 361a4d96a55SJed Brown @*/ 362a4d96a55SJed Brown PetscErrorCode ISLocalToGlobalMappingCreateSF(PetscSF sf,PetscInt start,ISLocalToGlobalMapping *mapping) 363a4d96a55SJed Brown { 364a4d96a55SJed Brown PetscErrorCode ierr; 365a4d96a55SJed Brown PetscInt i,maxlocal,nroots,nleaves,*globals,*ltog; 366a4d96a55SJed Brown const PetscInt *ilocal; 367a4d96a55SJed Brown MPI_Comm comm; 368a4d96a55SJed Brown 369a4d96a55SJed Brown PetscFunctionBegin; 370a4d96a55SJed Brown PetscValidHeaderSpecific(sf,PETSCSF_CLASSID,1); 371a4d96a55SJed Brown PetscValidPointer(mapping,3); 372a4d96a55SJed Brown 373a4d96a55SJed Brown ierr = PetscObjectGetComm((PetscObject)sf,&comm);CHKERRQ(ierr); 3740298fd71SBarry Smith ierr = PetscSFGetGraph(sf,&nroots,&nleaves,&ilocal,NULL);CHKERRQ(ierr); 375f6e5521dSKarl Rupp if (ilocal) { 376f6e5521dSKarl Rupp for (i=0,maxlocal=0; i<nleaves; i++) maxlocal = PetscMax(maxlocal,ilocal[i]+1); 377f6e5521dSKarl Rupp } 378a4d96a55SJed Brown else maxlocal = nleaves; 379785e854fSJed Brown ierr = PetscMalloc1(nroots,&globals);CHKERRQ(ierr); 380785e854fSJed Brown ierr = PetscMalloc1(maxlocal,<og);CHKERRQ(ierr); 381a4d96a55SJed Brown for (i=0; i<nroots; i++) globals[i] = start + i; 382a4d96a55SJed Brown for (i=0; i<maxlocal; i++) ltog[i] = -1; 383a4d96a55SJed Brown ierr = PetscSFBcastBegin(sf,MPIU_INT,globals,ltog);CHKERRQ(ierr); 384a4d96a55SJed Brown ierr = PetscSFBcastEnd(sf,MPIU_INT,globals,ltog);CHKERRQ(ierr); 385f0413b6fSBarry Smith ierr = ISLocalToGlobalMappingCreate(comm,1,maxlocal,ltog,PETSC_OWN_POINTER,mapping);CHKERRQ(ierr); 386a4d96a55SJed Brown ierr = PetscFree(globals);CHKERRQ(ierr); 387a4d96a55SJed Brown PetscFunctionReturn(0); 388a4d96a55SJed Brown } 389b46b645bSBarry Smith 39063fa5c83Sstefano_zampini /*@ 39163fa5c83Sstefano_zampini ISLocalToGlobalMappingSetBlockSize - Sets the blocksize of the mapping 39263fa5c83Sstefano_zampini 39363fa5c83Sstefano_zampini Not collective 39463fa5c83Sstefano_zampini 39563fa5c83Sstefano_zampini Input Parameters: 39663fa5c83Sstefano_zampini . mapping - mapping data structure 39763fa5c83Sstefano_zampini . bs - the blocksize 39863fa5c83Sstefano_zampini 39963fa5c83Sstefano_zampini Level: advanced 40063fa5c83Sstefano_zampini 40163fa5c83Sstefano_zampini Concepts: mapping^local to global 40263fa5c83Sstefano_zampini 40363fa5c83Sstefano_zampini .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS() 40463fa5c83Sstefano_zampini @*/ 40563fa5c83Sstefano_zampini PetscErrorCode ISLocalToGlobalMappingSetBlockSize(ISLocalToGlobalMapping mapping,PetscInt bs) 40663fa5c83Sstefano_zampini { 407a59f3c4dSstefano_zampini PetscInt *nid; 408a59f3c4dSstefano_zampini const PetscInt *oid; 409a59f3c4dSstefano_zampini PetscInt i,cn,on,obs,nn; 41063fa5c83Sstefano_zampini PetscErrorCode ierr; 41163fa5c83Sstefano_zampini 41263fa5c83Sstefano_zampini PetscFunctionBegin; 41363fa5c83Sstefano_zampini PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1); 41463fa5c83Sstefano_zampini if (bs < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid block size %D",bs); 41563fa5c83Sstefano_zampini if (bs == mapping->bs) PetscFunctionReturn(0); 41663fa5c83Sstefano_zampini on = mapping->n; 41763fa5c83Sstefano_zampini obs = mapping->bs; 41863fa5c83Sstefano_zampini oid = mapping->indices; 41963fa5c83Sstefano_zampini nn = (on*obs)/bs; 42063fa5c83Sstefano_zampini if ((on*obs)%bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Block size %D is inconsistent with block size %D and number of block indices %D",bs,obs,on); 421a59f3c4dSstefano_zampini 42263fa5c83Sstefano_zampini ierr = PetscMalloc1(nn,&nid);CHKERRQ(ierr); 423a59f3c4dSstefano_zampini ierr = ISLocalToGlobalMappingGetIndices(mapping,&oid);CHKERRQ(ierr); 424a59f3c4dSstefano_zampini for (i=0;i<nn;i++) { 425a59f3c4dSstefano_zampini PetscInt j; 426a59f3c4dSstefano_zampini for (j=0,cn=0;j<bs-1;j++) { 427a59f3c4dSstefano_zampini if (oid[i*bs+j] < 0) { cn++; continue; } 428a59f3c4dSstefano_zampini if (oid[i*bs+j] != oid[i*bs+j+1]-1) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Block sizes %D and %D are incompatible with the block indices: non consecutive indices %D %D",bs,obs,oid[i*bs+j],oid[i*bs+j+1]); 429a59f3c4dSstefano_zampini } 430a59f3c4dSstefano_zampini if (oid[i*bs+j] < 0) cn++; 4318b7cb0e6Sstefano_zampini if (cn) { 432a59f3c4dSstefano_zampini if (cn != bs) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Block sizes %D and %D are incompatible with the block indices: invalid number of negative entries in block %D",bs,obs,cn); 433a59f3c4dSstefano_zampini nid[i] = -1; 4348b7cb0e6Sstefano_zampini } else { 435a59f3c4dSstefano_zampini nid[i] = oid[i*bs]/bs; 43663fa5c83Sstefano_zampini } 43763fa5c83Sstefano_zampini } 438a59f3c4dSstefano_zampini ierr = ISLocalToGlobalMappingRestoreIndices(mapping,&oid);CHKERRQ(ierr); 439a59f3c4dSstefano_zampini 44063fa5c83Sstefano_zampini mapping->n = nn; 44163fa5c83Sstefano_zampini mapping->bs = bs; 44263fa5c83Sstefano_zampini ierr = PetscFree(mapping->indices);CHKERRQ(ierr); 44363fa5c83Sstefano_zampini mapping->indices = nid; 444c9345713Sstefano_zampini mapping->globalstart = 0; 445c9345713Sstefano_zampini mapping->globalend = 0; 4461bd0b88eSStefano Zampini 4471bd0b88eSStefano Zampini /* reset the cached information */ 4481bd0b88eSStefano Zampini ierr = PetscFree(mapping->info_procs);CHKERRQ(ierr); 4491bd0b88eSStefano Zampini ierr = PetscFree(mapping->info_numprocs);CHKERRQ(ierr); 4501bd0b88eSStefano Zampini if (mapping->info_indices) { 4511bd0b88eSStefano Zampini PetscInt i; 4521bd0b88eSStefano Zampini 4531bd0b88eSStefano Zampini ierr = PetscFree((mapping->info_indices)[0]);CHKERRQ(ierr); 4541bd0b88eSStefano Zampini for (i=1; i<mapping->info_nproc; i++) { 4551bd0b88eSStefano Zampini ierr = PetscFree(mapping->info_indices[i]);CHKERRQ(ierr); 4561bd0b88eSStefano Zampini } 4571bd0b88eSStefano Zampini ierr = PetscFree(mapping->info_indices);CHKERRQ(ierr); 4581bd0b88eSStefano Zampini } 4591bd0b88eSStefano Zampini mapping->info_cached = PETSC_FALSE; 4601bd0b88eSStefano Zampini 461413f72f0SBarry Smith if (mapping->ops->destroy) { 462413f72f0SBarry Smith ierr = (*mapping->ops->destroy)(mapping);CHKERRQ(ierr); 463413f72f0SBarry Smith } 46463fa5c83Sstefano_zampini PetscFunctionReturn(0); 46563fa5c83Sstefano_zampini } 46663fa5c83Sstefano_zampini 46745b6f7e9SBarry Smith /*@ 46845b6f7e9SBarry Smith ISLocalToGlobalMappingGetBlockSize - Gets the blocksize of the mapping 46945b6f7e9SBarry Smith ordering and a global parallel ordering. 47045b6f7e9SBarry Smith 47145b6f7e9SBarry Smith Not Collective 47245b6f7e9SBarry Smith 47345b6f7e9SBarry Smith Input Parameters: 47445b6f7e9SBarry Smith . mapping - mapping data structure 47545b6f7e9SBarry Smith 47645b6f7e9SBarry Smith Output Parameter: 47745b6f7e9SBarry Smith . bs - the blocksize 47845b6f7e9SBarry Smith 47945b6f7e9SBarry Smith Level: advanced 48045b6f7e9SBarry Smith 48145b6f7e9SBarry Smith Concepts: mapping^local to global 48245b6f7e9SBarry Smith 48345b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS() 48445b6f7e9SBarry Smith @*/ 48545b6f7e9SBarry Smith PetscErrorCode ISLocalToGlobalMappingGetBlockSize(ISLocalToGlobalMapping mapping,PetscInt *bs) 48645b6f7e9SBarry Smith { 48745b6f7e9SBarry Smith PetscFunctionBegin; 488cbc1caf0SMatthew G. Knepley PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1); 48945b6f7e9SBarry Smith *bs = mapping->bs; 49045b6f7e9SBarry Smith PetscFunctionReturn(0); 49145b6f7e9SBarry Smith } 49245b6f7e9SBarry Smith 493ba5bb76aSSatish Balay /*@ 49490f02eecSBarry Smith ISLocalToGlobalMappingCreate - Creates a mapping between a local (0 to n) 49590f02eecSBarry Smith ordering and a global parallel ordering. 4962362add9SBarry Smith 49789d82c54SBarry Smith Not Collective, but communicator may have more than one process 498b9cd556bSLois Curfman McInnes 4992362add9SBarry Smith Input Parameters: 50089d82c54SBarry Smith + comm - MPI communicator 501f0413b6fSBarry Smith . bs - the block size 50228bc9809SBarry Smith . n - the number of local elements divided by the block size, or equivalently the number of block indices 50328bc9809SBarry 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 504d5ad8652SBarry Smith - mode - see PetscCopyMode 5052362add9SBarry Smith 506a997ad1aSLois Curfman McInnes Output Parameter: 50790f02eecSBarry Smith . mapping - new mapping data structure 5082362add9SBarry Smith 50995452b02SPatrick Sanan Notes: 51095452b02SPatrick Sanan There is one integer value in indices per block and it represents the actual indices bs*idx + j, where j=0,..,bs-1 511413f72f0SBarry Smith 5129a7b7924SJed Brown For "small" problems when using ISGlobalToLocalMappingApply() and ISGlobalToLocalMappingApplyBlock(), the ISLocalToGlobalMappingType of ISLOCALTOGLOBALMAPPINGBASIC will be used; 513413f72f0SBarry 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. 514413f72f0SBarry Smith Use ISLocalToGlobalMappingSetType() or call ISLocalToGlobalMappingSetFromOptions() with the option -islocaltoglobalmapping_type <basic,hash> to control which is used. 515413f72f0SBarry Smith 516a997ad1aSLois Curfman McInnes Level: advanced 517a997ad1aSLois Curfman McInnes 518273d9f13SBarry Smith Concepts: mapping^local to global 5192362add9SBarry Smith 520413f72f0SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingSetFromOptions(), ISLOCALTOGLOBALMAPPINGBASIC, ISLOCALTOGLOBALMAPPINGHASH 521413f72f0SBarry Smith ISLocalToGlobalMappingSetType(), ISLocalToGlobalMappingType 5222362add9SBarry Smith @*/ 52360c7cefcSBarry Smith PetscErrorCode ISLocalToGlobalMappingCreate(MPI_Comm comm,PetscInt bs,PetscInt n,const PetscInt indices[],PetscCopyMode mode,ISLocalToGlobalMapping *mapping) 5242362add9SBarry Smith { 5256849ba73SBarry Smith PetscErrorCode ierr; 52632dcc486SBarry Smith PetscInt *in; 527b46b645bSBarry Smith 528b46b645bSBarry Smith PetscFunctionBegin; 52973911063SBarry Smith if (n) PetscValidIntPointer(indices,3); 5304482741eSBarry Smith PetscValidPointer(mapping,4); 531b46b645bSBarry Smith 5320298fd71SBarry Smith *mapping = NULL; 533607a6623SBarry Smith ierr = ISInitializePackage();CHKERRQ(ierr); 5342362add9SBarry Smith 53573107ff1SLisandro Dalcin ierr = PetscHeaderCreate(*mapping,IS_LTOGM_CLASSID,"ISLocalToGlobalMapping","Local to global mapping","IS", 53660c7cefcSBarry Smith comm,ISLocalToGlobalMappingDestroy,ISLocalToGlobalMappingView);CHKERRQ(ierr); 537d4bb536fSBarry Smith (*mapping)->n = n; 538f0413b6fSBarry Smith (*mapping)->bs = bs; 539268a049cSStefano Zampini (*mapping)->info_cached = PETSC_FALSE; 540268a049cSStefano Zampini (*mapping)->info_free = PETSC_FALSE; 541268a049cSStefano Zampini (*mapping)->info_procs = NULL; 542268a049cSStefano Zampini (*mapping)->info_numprocs = NULL; 543268a049cSStefano Zampini (*mapping)->info_indices = NULL; 5441bd0b88eSStefano Zampini (*mapping)->info_nodec = NULL; 5451bd0b88eSStefano Zampini (*mapping)->info_nodei = NULL; 546413f72f0SBarry Smith 547413f72f0SBarry Smith (*mapping)->ops->globaltolocalmappingapply = NULL; 548413f72f0SBarry Smith (*mapping)->ops->globaltolocalmappingapplyblock = NULL; 549413f72f0SBarry Smith (*mapping)->ops->destroy = NULL; 550d5ad8652SBarry Smith if (mode == PETSC_COPY_VALUES) { 551785e854fSJed Brown ierr = PetscMalloc1(n,&in);CHKERRQ(ierr); 552d5ad8652SBarry Smith ierr = PetscMemcpy(in,indices,n*sizeof(PetscInt));CHKERRQ(ierr); 553d5ad8652SBarry Smith (*mapping)->indices = in; 5546389a1a1SBarry Smith ierr = PetscLogObjectMemory((PetscObject)*mapping,n*sizeof(PetscInt));CHKERRQ(ierr); 5556389a1a1SBarry Smith } else if (mode == PETSC_OWN_POINTER) { 5566389a1a1SBarry Smith (*mapping)->indices = (PetscInt*)indices; 5576389a1a1SBarry Smith ierr = PetscLogObjectMemory((PetscObject)*mapping,n*sizeof(PetscInt));CHKERRQ(ierr); 5586389a1a1SBarry Smith } 55960c7cefcSBarry Smith else SETERRQ(comm,PETSC_ERR_SUP,"Cannot currently use PETSC_USE_POINTER"); 5603a40ed3dSBarry Smith PetscFunctionReturn(0); 5612362add9SBarry Smith } 5622362add9SBarry Smith 563413f72f0SBarry Smith PetscFunctionList ISLocalToGlobalMappingList = NULL; 564413f72f0SBarry Smith 565413f72f0SBarry Smith 56690f02eecSBarry Smith /*@ 5677e99dc12SLawrence Mitchell ISLocalToGlobalMappingSetFromOptions - Set mapping options from the options database. 5687e99dc12SLawrence Mitchell 5697e99dc12SLawrence Mitchell Not collective 5707e99dc12SLawrence Mitchell 5717e99dc12SLawrence Mitchell Input Parameters: 5727e99dc12SLawrence Mitchell . mapping - mapping data structure 5737e99dc12SLawrence Mitchell 5747e99dc12SLawrence Mitchell Level: advanced 5757e99dc12SLawrence Mitchell 5767e99dc12SLawrence Mitchell @*/ 5777e99dc12SLawrence Mitchell PetscErrorCode ISLocalToGlobalMappingSetFromOptions(ISLocalToGlobalMapping mapping) 5787e99dc12SLawrence Mitchell { 5797e99dc12SLawrence Mitchell PetscErrorCode ierr; 580413f72f0SBarry Smith char type[256]; 581413f72f0SBarry Smith ISLocalToGlobalMappingType defaulttype = "Not set"; 5827e99dc12SLawrence Mitchell PetscBool flg; 5837e99dc12SLawrence Mitchell 5847e99dc12SLawrence Mitchell PetscFunctionBegin; 5857e99dc12SLawrence Mitchell PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1); 586413f72f0SBarry Smith ierr = ISLocalToGlobalMappingRegisterAll();CHKERRQ(ierr); 5877e99dc12SLawrence Mitchell ierr = PetscObjectOptionsBegin((PetscObject)mapping);CHKERRQ(ierr); 588413f72f0SBarry Smith ierr = PetscOptionsFList("-islocaltoglobalmapping_type","ISLocalToGlobalMapping method","ISLocalToGlobalMappingSetType",ISLocalToGlobalMappingList,(char*)(((PetscObject)mapping)->type_name) ? ((PetscObject)mapping)->type_name : defaulttype,type,256,&flg);CHKERRQ(ierr); 589413f72f0SBarry Smith if (flg) { 590413f72f0SBarry Smith ierr = ISLocalToGlobalMappingSetType(mapping,type);CHKERRQ(ierr); 591413f72f0SBarry Smith } 5927e99dc12SLawrence Mitchell ierr = PetscOptionsEnd();CHKERRQ(ierr); 5937e99dc12SLawrence Mitchell PetscFunctionReturn(0); 5947e99dc12SLawrence Mitchell } 5957e99dc12SLawrence Mitchell 5967e99dc12SLawrence Mitchell /*@ 59790f02eecSBarry Smith ISLocalToGlobalMappingDestroy - Destroys a mapping between a local (0 to n) 59890f02eecSBarry Smith ordering and a global parallel ordering. 59990f02eecSBarry Smith 6000f5bd95cSBarry Smith Note Collective 601b9cd556bSLois Curfman McInnes 60290f02eecSBarry Smith Input Parameters: 60390f02eecSBarry Smith . mapping - mapping data structure 60490f02eecSBarry Smith 605a997ad1aSLois Curfman McInnes Level: advanced 606a997ad1aSLois Curfman McInnes 6073acfe500SLois Curfman McInnes .seealso: ISLocalToGlobalMappingCreate() 60890f02eecSBarry Smith @*/ 6096bf464f9SBarry Smith PetscErrorCode ISLocalToGlobalMappingDestroy(ISLocalToGlobalMapping *mapping) 61090f02eecSBarry Smith { 611dfbe8321SBarry Smith PetscErrorCode ierr; 6125fd66863SKarl Rupp 6133a40ed3dSBarry Smith PetscFunctionBegin; 6146bf464f9SBarry Smith if (!*mapping) PetscFunctionReturn(0); 6156bf464f9SBarry Smith PetscValidHeaderSpecific((*mapping),IS_LTOGM_CLASSID,1); 616997056adSBarry Smith if (--((PetscObject)(*mapping))->refct > 0) {*mapping = 0;PetscFunctionReturn(0);} 6176bf464f9SBarry Smith ierr = PetscFree((*mapping)->indices);CHKERRQ(ierr); 618268a049cSStefano Zampini ierr = PetscFree((*mapping)->info_procs);CHKERRQ(ierr); 619268a049cSStefano Zampini ierr = PetscFree((*mapping)->info_numprocs);CHKERRQ(ierr); 620268a049cSStefano Zampini if ((*mapping)->info_indices) { 621268a049cSStefano Zampini PetscInt i; 622268a049cSStefano Zampini 623268a049cSStefano Zampini ierr = PetscFree(((*mapping)->info_indices)[0]);CHKERRQ(ierr); 624268a049cSStefano Zampini for (i=1; i<(*mapping)->info_nproc; i++) { 625268a049cSStefano Zampini ierr = PetscFree(((*mapping)->info_indices)[i]);CHKERRQ(ierr); 626268a049cSStefano Zampini } 627268a049cSStefano Zampini ierr = PetscFree((*mapping)->info_indices);CHKERRQ(ierr); 628268a049cSStefano Zampini } 6291bd0b88eSStefano Zampini if ((*mapping)->info_nodei) { 6301bd0b88eSStefano Zampini ierr = PetscFree(((*mapping)->info_nodei)[0]);CHKERRQ(ierr); 6311bd0b88eSStefano Zampini } 6321bd0b88eSStefano Zampini ierr = PetscFree((*mapping)->info_nodei);CHKERRQ(ierr); 6331bd0b88eSStefano Zampini ierr = PetscFree((*mapping)->info_nodec);CHKERRQ(ierr); 634413f72f0SBarry Smith if ((*mapping)->ops->destroy) { 635413f72f0SBarry Smith ierr = (*(*mapping)->ops->destroy)(*mapping);CHKERRQ(ierr); 636413f72f0SBarry Smith } 637d38fa0fbSBarry Smith ierr = PetscHeaderDestroy(mapping);CHKERRQ(ierr); 638992144d0SBarry Smith *mapping = 0; 6393a40ed3dSBarry Smith PetscFunctionReturn(0); 64090f02eecSBarry Smith } 64190f02eecSBarry Smith 64290f02eecSBarry Smith /*@ 6433acfe500SLois Curfman McInnes ISLocalToGlobalMappingApplyIS - Creates from an IS in the local numbering 6443acfe500SLois Curfman McInnes a new index set using the global numbering defined in an ISLocalToGlobalMapping 6453acfe500SLois Curfman McInnes context. 64690f02eecSBarry Smith 647b9cd556bSLois Curfman McInnes Not collective 648b9cd556bSLois Curfman McInnes 64990f02eecSBarry Smith Input Parameters: 650b9cd556bSLois Curfman McInnes + mapping - mapping between local and global numbering 651b9cd556bSLois Curfman McInnes - is - index set in local numbering 65290f02eecSBarry Smith 65390f02eecSBarry Smith Output Parameters: 65490f02eecSBarry Smith . newis - index set in global numbering 65590f02eecSBarry Smith 656a997ad1aSLois Curfman McInnes Level: advanced 657a997ad1aSLois Curfman McInnes 658273d9f13SBarry Smith Concepts: mapping^local to global 6593acfe500SLois Curfman McInnes 66090f02eecSBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingCreate(), 661d4bb536fSBarry Smith ISLocalToGlobalMappingDestroy(), ISGlobalToLocalMappingApply() 66290f02eecSBarry Smith @*/ 6637087cfbeSBarry Smith PetscErrorCode ISLocalToGlobalMappingApplyIS(ISLocalToGlobalMapping mapping,IS is,IS *newis) 66490f02eecSBarry Smith { 6656849ba73SBarry Smith PetscErrorCode ierr; 666e24637baSBarry Smith PetscInt n,*idxout; 6675d0c19d7SBarry Smith const PetscInt *idxin; 6683a40ed3dSBarry Smith 6693a40ed3dSBarry Smith PetscFunctionBegin; 6700700a824SBarry Smith PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1); 6710700a824SBarry Smith PetscValidHeaderSpecific(is,IS_CLASSID,2); 6724482741eSBarry Smith PetscValidPointer(newis,3); 67390f02eecSBarry Smith 6743b9aefa3SBarry Smith ierr = ISGetLocalSize(is,&n);CHKERRQ(ierr); 67590f02eecSBarry Smith ierr = ISGetIndices(is,&idxin);CHKERRQ(ierr); 676785e854fSJed Brown ierr = PetscMalloc1(n,&idxout);CHKERRQ(ierr); 677e24637baSBarry Smith ierr = ISLocalToGlobalMappingApply(mapping,n,idxin,idxout);CHKERRQ(ierr); 6783b9aefa3SBarry Smith ierr = ISRestoreIndices(is,&idxin);CHKERRQ(ierr); 679543f3098SMatthew G. Knepley ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is),n,idxout,PETSC_OWN_POINTER,newis);CHKERRQ(ierr); 6803a40ed3dSBarry Smith PetscFunctionReturn(0); 68190f02eecSBarry Smith } 68290f02eecSBarry Smith 683b89cb25eSSatish Balay /*@ 6843acfe500SLois Curfman McInnes ISLocalToGlobalMappingApply - Takes a list of integers in a local numbering 6853acfe500SLois Curfman McInnes and converts them to the global numbering. 68690f02eecSBarry Smith 687b9cd556bSLois Curfman McInnes Not collective 688b9cd556bSLois Curfman McInnes 689bb25748dSBarry Smith Input Parameters: 690b9cd556bSLois Curfman McInnes + mapping - the local to global mapping context 691bb25748dSBarry Smith . N - number of integers 692b9cd556bSLois Curfman McInnes - in - input indices in local numbering 693bb25748dSBarry Smith 694bb25748dSBarry Smith Output Parameter: 695bb25748dSBarry Smith . out - indices in global numbering 696bb25748dSBarry Smith 697b9cd556bSLois Curfman McInnes Notes: 698b9cd556bSLois Curfman McInnes The in and out array parameters may be identical. 699d4bb536fSBarry Smith 700a997ad1aSLois Curfman McInnes Level: advanced 701a997ad1aSLois Curfman McInnes 70245b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingApplyBlock(), ISLocalToGlobalMappingCreate(),ISLocalToGlobalMappingDestroy(), 7030752156aSBarry Smith ISLocalToGlobalMappingApplyIS(),AOCreateBasic(),AOApplicationToPetsc(), 704d4bb536fSBarry Smith AOPetscToApplication(), ISGlobalToLocalMappingApply() 705bb25748dSBarry Smith 706273d9f13SBarry Smith Concepts: mapping^local to global 707afcb2eb5SJed Brown @*/ 708afcb2eb5SJed Brown PetscErrorCode ISLocalToGlobalMappingApply(ISLocalToGlobalMapping mapping,PetscInt N,const PetscInt in[],PetscInt out[]) 709afcb2eb5SJed Brown { 710cbc1caf0SMatthew G. Knepley PetscInt i,bs,Nmax; 71145b6f7e9SBarry Smith 71245b6f7e9SBarry Smith PetscFunctionBegin; 713cbc1caf0SMatthew G. Knepley PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1); 714cbc1caf0SMatthew G. Knepley bs = mapping->bs; 715cbc1caf0SMatthew G. Knepley Nmax = bs*mapping->n; 71645b6f7e9SBarry Smith if (bs == 1) { 717cbc1caf0SMatthew G. Knepley const PetscInt *idx = mapping->indices; 71845b6f7e9SBarry Smith for (i=0; i<N; i++) { 71945b6f7e9SBarry Smith if (in[i] < 0) { 72045b6f7e9SBarry Smith out[i] = in[i]; 72145b6f7e9SBarry Smith continue; 72245b6f7e9SBarry Smith } 723e24637baSBarry Smith if (in[i] >= Nmax) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local index %D too large %D (max) at %D",in[i],Nmax-1,i); 72445b6f7e9SBarry Smith out[i] = idx[in[i]]; 72545b6f7e9SBarry Smith } 72645b6f7e9SBarry Smith } else { 727cbc1caf0SMatthew G. Knepley const PetscInt *idx = mapping->indices; 72845b6f7e9SBarry Smith for (i=0; i<N; i++) { 72945b6f7e9SBarry Smith if (in[i] < 0) { 73045b6f7e9SBarry Smith out[i] = in[i]; 73145b6f7e9SBarry Smith continue; 73245b6f7e9SBarry Smith } 733e24637baSBarry Smith if (in[i] >= Nmax) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local index %D too large %D (max) at %D",in[i],Nmax-1,i); 73445b6f7e9SBarry Smith out[i] = idx[in[i]/bs]*bs + (in[i] % bs); 73545b6f7e9SBarry Smith } 73645b6f7e9SBarry Smith } 73745b6f7e9SBarry Smith PetscFunctionReturn(0); 73845b6f7e9SBarry Smith } 73945b6f7e9SBarry Smith 74045b6f7e9SBarry Smith /*@ 7416006e8d2SBarry Smith ISLocalToGlobalMappingApplyBlock - Takes a list of integers in a local block numbering and converts them to the global block numbering 74245b6f7e9SBarry Smith 74345b6f7e9SBarry Smith Not collective 74445b6f7e9SBarry Smith 74545b6f7e9SBarry Smith Input Parameters: 74645b6f7e9SBarry Smith + mapping - the local to global mapping context 74745b6f7e9SBarry Smith . N - number of integers 7486006e8d2SBarry Smith - in - input indices in local block numbering 74945b6f7e9SBarry Smith 75045b6f7e9SBarry Smith Output Parameter: 7516006e8d2SBarry Smith . out - indices in global block numbering 75245b6f7e9SBarry Smith 75345b6f7e9SBarry Smith Notes: 75445b6f7e9SBarry Smith The in and out array parameters may be identical. 75545b6f7e9SBarry Smith 7566006e8d2SBarry Smith Example: 7576006e8d2SBarry 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 7586006e8d2SBarry Smith (the first block) would produce 0 and the mapping applied to 1 (the second block) would produce 3. 7596006e8d2SBarry Smith 76045b6f7e9SBarry Smith Level: advanced 76145b6f7e9SBarry Smith 76245b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingCreate(),ISLocalToGlobalMappingDestroy(), 76345b6f7e9SBarry Smith ISLocalToGlobalMappingApplyIS(),AOCreateBasic(),AOApplicationToPetsc(), 76445b6f7e9SBarry Smith AOPetscToApplication(), ISGlobalToLocalMappingApply() 76545b6f7e9SBarry Smith 76645b6f7e9SBarry Smith Concepts: mapping^local to global 76745b6f7e9SBarry Smith @*/ 76845b6f7e9SBarry Smith PetscErrorCode ISLocalToGlobalMappingApplyBlock(ISLocalToGlobalMapping mapping,PetscInt N,const PetscInt in[],PetscInt out[]) 76945b6f7e9SBarry Smith { 770cbc1caf0SMatthew G. Knepley 771cbc1caf0SMatthew G. Knepley PetscFunctionBegin; 772cbc1caf0SMatthew G. Knepley PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1); 773cbc1caf0SMatthew G. Knepley { 774afcb2eb5SJed Brown PetscInt i,Nmax = mapping->n; 775afcb2eb5SJed Brown const PetscInt *idx = mapping->indices; 776d4bb536fSBarry Smith 777afcb2eb5SJed Brown for (i=0; i<N; i++) { 778afcb2eb5SJed Brown if (in[i] < 0) { 779afcb2eb5SJed Brown out[i] = in[i]; 780afcb2eb5SJed Brown continue; 781afcb2eb5SJed Brown } 782e24637baSBarry Smith if (in[i] >= Nmax) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Local block index %D too large %D (max) at %D",in[i],Nmax-1,i); 783afcb2eb5SJed Brown out[i] = idx[in[i]]; 784afcb2eb5SJed Brown } 785cbc1caf0SMatthew G. Knepley } 786afcb2eb5SJed Brown PetscFunctionReturn(0); 787afcb2eb5SJed Brown } 788d4bb536fSBarry Smith 7897e99dc12SLawrence Mitchell /*@ 790a997ad1aSLois Curfman McInnes ISGlobalToLocalMappingApply - Provides the local numbering for a list of integers 791a997ad1aSLois Curfman McInnes specified with a global numbering. 792d4bb536fSBarry Smith 793b9cd556bSLois Curfman McInnes Not collective 794b9cd556bSLois Curfman McInnes 795d4bb536fSBarry Smith Input Parameters: 796b9cd556bSLois Curfman McInnes + mapping - mapping between local and global numbering 797d4bb536fSBarry Smith . type - IS_GTOLM_MASK - replaces global indices with no local value with -1 798d4bb536fSBarry Smith IS_GTOLM_DROP - drops the indices with no local value from the output list 799d4bb536fSBarry Smith . n - number of global indices to map 800b9cd556bSLois Curfman McInnes - idx - global indices to map 801d4bb536fSBarry Smith 802d4bb536fSBarry Smith Output Parameters: 803b9cd556bSLois Curfman McInnes + nout - number of indices in output array (if type == IS_GTOLM_MASK then nout = n) 804b9cd556bSLois Curfman McInnes - idxout - local index of each global index, one must pass in an array long enough 805e182c471SBarry Smith to hold all the indices. You can call ISGlobalToLocalMappingApply() with 8060298fd71SBarry Smith idxout == NULL to determine the required length (returned in nout) 807e182c471SBarry Smith and then allocate the required space and call ISGlobalToLocalMappingApply() 808e182c471SBarry Smith a second time to set the values. 809d4bb536fSBarry Smith 810b9cd556bSLois Curfman McInnes Notes: 8110298fd71SBarry Smith Either nout or idxout may be NULL. idx and idxout may be identical. 812d4bb536fSBarry Smith 8139a7b7924SJed Brown For "small" problems when using ISGlobalToLocalMappingApply() and ISGlobalToLocalMappingApplyBlock(), the ISLocalToGlobalMappingType of ISLOCALTOGLOBALMAPPINGBASIC will be used; 814413f72f0SBarry 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. 815413f72f0SBarry Smith Use ISLocalToGlobalMappingSetType() or call ISLocalToGlobalMappingSetFromOptions() with the option -islocaltoglobalmapping_type <basic,hash> to control which is used. 8160f5bd95cSBarry Smith 817a997ad1aSLois Curfman McInnes Level: advanced 818a997ad1aSLois Curfman McInnes 81932fd6b96SBarry Smith Developer Note: The manual page states that idx and idxout may be identical but the calling 82032fd6b96SBarry Smith sequence declares idx as const so it cannot be the same as idxout. 82132fd6b96SBarry Smith 822273d9f13SBarry Smith Concepts: mapping^global to local 823d4bb536fSBarry Smith 8249d90f715SBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISGlobalToLocalMappingApplyBlock(), ISLocalToGlobalMappingCreate(), 825413f72f0SBarry Smith ISLocalToGlobalMappingDestroy() 826d4bb536fSBarry Smith @*/ 827413f72f0SBarry Smith PetscErrorCode ISGlobalToLocalMappingApply(ISLocalToGlobalMapping mapping,ISGlobalToLocalMappingMode type,PetscInt n,const PetscInt idx[],PetscInt *nout,PetscInt idxout[]) 828d4bb536fSBarry Smith { 8299d90f715SBarry Smith PetscErrorCode ierr; 8309d90f715SBarry Smith 8319d90f715SBarry Smith PetscFunctionBegin; 8329d90f715SBarry Smith PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1); 833413f72f0SBarry Smith if (!mapping->data) { 834413f72f0SBarry Smith ierr = ISGlobalToLocalMappingSetUp(mapping);CHKERRQ(ierr); 8359d90f715SBarry Smith } 836413f72f0SBarry Smith ierr = (*mapping->ops->globaltolocalmappingapply)(mapping,type,n,idx,nout,idxout);CHKERRQ(ierr); 8379d90f715SBarry Smith PetscFunctionReturn(0); 8389d90f715SBarry Smith } 8399d90f715SBarry Smith 840d4fe737eSStefano Zampini /*@ 841d4fe737eSStefano Zampini ISGlobalToLocalMappingApplyIS - Creates from an IS in the global numbering 842d4fe737eSStefano Zampini a new index set using the local numbering defined in an ISLocalToGlobalMapping 843d4fe737eSStefano Zampini context. 844d4fe737eSStefano Zampini 845d4fe737eSStefano Zampini Not collective 846d4fe737eSStefano Zampini 847d4fe737eSStefano Zampini Input Parameters: 848d4fe737eSStefano Zampini + mapping - mapping between local and global numbering 8492785b321SStefano Zampini . type - IS_GTOLM_MASK - replaces global indices with no local value with -1 8502785b321SStefano Zampini IS_GTOLM_DROP - drops the indices with no local value from the output list 851d4fe737eSStefano Zampini - is - index set in global numbering 852d4fe737eSStefano Zampini 853d4fe737eSStefano Zampini Output Parameters: 854d4fe737eSStefano Zampini . newis - index set in local numbering 855d4fe737eSStefano Zampini 856d4fe737eSStefano Zampini Level: advanced 857d4fe737eSStefano Zampini 858d4fe737eSStefano Zampini Concepts: mapping^local to global 859d4fe737eSStefano Zampini 860d4fe737eSStefano Zampini .seealso: ISGlobalToLocalMappingApply(), ISLocalToGlobalMappingCreate(), 861d4fe737eSStefano Zampini ISLocalToGlobalMappingDestroy() 862d4fe737eSStefano Zampini @*/ 863413f72f0SBarry Smith PetscErrorCode ISGlobalToLocalMappingApplyIS(ISLocalToGlobalMapping mapping,ISGlobalToLocalMappingMode type,IS is,IS *newis) 864d4fe737eSStefano Zampini { 865d4fe737eSStefano Zampini PetscErrorCode ierr; 866d4fe737eSStefano Zampini PetscInt n,nout,*idxout; 867d4fe737eSStefano Zampini const PetscInt *idxin; 868d4fe737eSStefano Zampini 869d4fe737eSStefano Zampini PetscFunctionBegin; 870d4fe737eSStefano Zampini PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1); 871d4fe737eSStefano Zampini PetscValidHeaderSpecific(is,IS_CLASSID,3); 872d4fe737eSStefano Zampini PetscValidPointer(newis,4); 873d4fe737eSStefano Zampini 874d4fe737eSStefano Zampini ierr = ISGetLocalSize(is,&n);CHKERRQ(ierr); 875d4fe737eSStefano Zampini ierr = ISGetIndices(is,&idxin);CHKERRQ(ierr); 876d4fe737eSStefano Zampini if (type == IS_GTOLM_MASK) { 877d4fe737eSStefano Zampini ierr = PetscMalloc1(n,&idxout);CHKERRQ(ierr); 878d4fe737eSStefano Zampini } else { 879d4fe737eSStefano Zampini ierr = ISGlobalToLocalMappingApply(mapping,type,n,idxin,&nout,NULL);CHKERRQ(ierr); 880d4fe737eSStefano Zampini ierr = PetscMalloc1(nout,&idxout);CHKERRQ(ierr); 881d4fe737eSStefano Zampini } 882d4fe737eSStefano Zampini ierr = ISGlobalToLocalMappingApply(mapping,type,n,idxin,&nout,idxout);CHKERRQ(ierr); 883d4fe737eSStefano Zampini ierr = ISRestoreIndices(is,&idxin);CHKERRQ(ierr); 884d4fe737eSStefano Zampini ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is),nout,idxout,PETSC_OWN_POINTER,newis);CHKERRQ(ierr); 885d4fe737eSStefano Zampini PetscFunctionReturn(0); 886d4fe737eSStefano Zampini } 887d4fe737eSStefano Zampini 8889d90f715SBarry Smith /*@ 8899d90f715SBarry Smith ISGlobalToLocalMappingApplyBlock - Provides the local block numbering for a list of integers 8909d90f715SBarry Smith specified with a block global numbering. 8919d90f715SBarry Smith 8929d90f715SBarry Smith Not collective 8939d90f715SBarry Smith 8949d90f715SBarry Smith Input Parameters: 8959d90f715SBarry Smith + mapping - mapping between local and global numbering 8969d90f715SBarry Smith . type - IS_GTOLM_MASK - replaces global indices with no local value with -1 8979d90f715SBarry Smith IS_GTOLM_DROP - drops the indices with no local value from the output list 8989d90f715SBarry Smith . n - number of global indices to map 8999d90f715SBarry Smith - idx - global indices to map 9009d90f715SBarry Smith 9019d90f715SBarry Smith Output Parameters: 9029d90f715SBarry Smith + nout - number of indices in output array (if type == IS_GTOLM_MASK then nout = n) 9039d90f715SBarry Smith - idxout - local index of each global index, one must pass in an array long enough 9049d90f715SBarry Smith to hold all the indices. You can call ISGlobalToLocalMappingApplyBlock() with 9059d90f715SBarry Smith idxout == NULL to determine the required length (returned in nout) 9069d90f715SBarry Smith and then allocate the required space and call ISGlobalToLocalMappingApplyBlock() 9079d90f715SBarry Smith a second time to set the values. 9089d90f715SBarry Smith 9099d90f715SBarry Smith Notes: 9109d90f715SBarry Smith Either nout or idxout may be NULL. idx and idxout may be identical. 9119d90f715SBarry Smith 9129a7b7924SJed Brown For "small" problems when using ISGlobalToLocalMappingApply() and ISGlobalToLocalMappingApplyBlock(), the ISLocalToGlobalMappingType of ISLOCALTOGLOBALMAPPINGBASIC will be used; 913413f72f0SBarry 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. 914413f72f0SBarry Smith Use ISLocalToGlobalMappingSetType() or call ISLocalToGlobalMappingSetFromOptions() with the option -islocaltoglobalmapping_type <basic,hash> to control which is used. 9159d90f715SBarry Smith 9169d90f715SBarry Smith Level: advanced 9179d90f715SBarry Smith 9189d90f715SBarry Smith Developer Note: The manual page states that idx and idxout may be identical but the calling 9199d90f715SBarry Smith sequence declares idx as const so it cannot be the same as idxout. 9209d90f715SBarry Smith 9219d90f715SBarry Smith Concepts: mapping^global to local 9229d90f715SBarry Smith 9239d90f715SBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISGlobalToLocalMappingApply(), ISLocalToGlobalMappingCreate(), 924413f72f0SBarry Smith ISLocalToGlobalMappingDestroy() 9259d90f715SBarry Smith @*/ 926413f72f0SBarry Smith PetscErrorCode ISGlobalToLocalMappingApplyBlock(ISLocalToGlobalMapping mapping,ISGlobalToLocalMappingMode type, 9279d90f715SBarry Smith PetscInt n,const PetscInt idx[],PetscInt *nout,PetscInt idxout[]) 9289d90f715SBarry Smith { 9296849ba73SBarry Smith PetscErrorCode ierr; 930d4bb536fSBarry Smith 9313a40ed3dSBarry Smith PetscFunctionBegin; 9320700a824SBarry Smith PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1); 933413f72f0SBarry Smith if (!mapping->data) { 934413f72f0SBarry Smith ierr = ISGlobalToLocalMappingSetUp(mapping);CHKERRQ(ierr); 935d4bb536fSBarry Smith } 936413f72f0SBarry Smith ierr = (*mapping->ops->globaltolocalmappingapplyblock)(mapping,type,n,idx,nout,idxout);CHKERRQ(ierr); 9373a40ed3dSBarry Smith PetscFunctionReturn(0); 938d4bb536fSBarry Smith } 93990f02eecSBarry Smith 940413f72f0SBarry Smith 94189d82c54SBarry Smith /*@C 9426a818285SBarry Smith ISLocalToGlobalMappingGetBlockInfo - Gets the neighbor information for each processor and 94389d82c54SBarry Smith each index shared by more than one processor 94489d82c54SBarry Smith 94589d82c54SBarry Smith Collective on ISLocalToGlobalMapping 94689d82c54SBarry Smith 94789d82c54SBarry Smith Input Parameters: 94889d82c54SBarry Smith . mapping - the mapping from local to global indexing 94989d82c54SBarry Smith 95089d82c54SBarry Smith Output Parameter: 95189d82c54SBarry Smith + nproc - number of processors that are connected to this one 95289d82c54SBarry Smith . proc - neighboring processors 95307b52d57SBarry Smith . numproc - number of indices for each subdomain (processor) 9543463a7baSJed Brown - indices - indices of nodes (in local numbering) shared with neighbors (sorted by global numbering) 95589d82c54SBarry Smith 95689d82c54SBarry Smith Level: advanced 95789d82c54SBarry Smith 958273d9f13SBarry Smith Concepts: mapping^local to global 95989d82c54SBarry Smith 9602cfcea29SBarry Smith Fortran Usage: 9612cfcea29SBarry Smith $ ISLocalToGlobalMpngGetInfoSize(ISLocalToGlobalMapping,PetscInt nproc,PetscInt numprocmax,ierr) followed by 9622cfcea29SBarry Smith $ ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt nproc, PetscInt procs[nproc],PetscInt numprocs[nproc], 9632cfcea29SBarry Smith PetscInt indices[nproc][numprocmax],ierr) 9642cfcea29SBarry Smith There is no ISLocalToGlobalMappingRestoreInfo() in Fortran. You must make sure that procs[], numprocs[] and 9652cfcea29SBarry Smith indices[][] are large enough arrays, either by allocating them dynamically or defining static ones large enough. 9662cfcea29SBarry Smith 9672cfcea29SBarry Smith 96807b52d57SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(), 96907b52d57SBarry Smith ISLocalToGlobalMappingRestoreInfo() 97089d82c54SBarry Smith @*/ 9716a818285SBarry Smith PetscErrorCode ISLocalToGlobalMappingGetBlockInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[]) 97289d82c54SBarry Smith { 9736849ba73SBarry Smith PetscErrorCode ierr; 974268a049cSStefano Zampini 975268a049cSStefano Zampini PetscFunctionBegin; 976268a049cSStefano Zampini PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1); 977268a049cSStefano Zampini if (mapping->info_cached) { 978268a049cSStefano Zampini *nproc = mapping->info_nproc; 979268a049cSStefano Zampini *procs = mapping->info_procs; 980268a049cSStefano Zampini *numprocs = mapping->info_numprocs; 981268a049cSStefano Zampini *indices = mapping->info_indices; 982268a049cSStefano Zampini } else { 983268a049cSStefano Zampini ierr = ISLocalToGlobalMappingGetBlockInfo_Private(mapping,nproc,procs,numprocs,indices);CHKERRQ(ierr); 984268a049cSStefano Zampini } 985268a049cSStefano Zampini PetscFunctionReturn(0); 986268a049cSStefano Zampini } 987268a049cSStefano Zampini 988268a049cSStefano Zampini static PetscErrorCode ISLocalToGlobalMappingGetBlockInfo_Private(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[]) 989268a049cSStefano Zampini { 990268a049cSStefano Zampini PetscErrorCode ierr; 99197f1f81fSBarry Smith PetscMPIInt size,rank,tag1,tag2,tag3,*len,*source,imdex; 99232dcc486SBarry Smith PetscInt i,n = mapping->n,Ng,ng,max = 0,*lindices = mapping->indices; 99332dcc486SBarry Smith PetscInt *nprocs,*owner,nsends,*sends,j,*starts,nmax,nrecvs,*recvs,proc; 99497f1f81fSBarry Smith PetscInt cnt,scale,*ownedsenders,*nownedsenders,rstart,nowned; 99532dcc486SBarry Smith PetscInt node,nownedm,nt,*sends2,nsends2,*starts2,*lens2,*dest,nrecvs2,*starts3,*recvs2,k,*bprocs,*tmp; 99632dcc486SBarry Smith PetscInt first_procs,first_numprocs,*first_indices; 99789d82c54SBarry Smith MPI_Request *recv_waits,*send_waits; 99830dcb7c9SBarry Smith MPI_Status recv_status,*send_status,*recv_statuses; 999ce94432eSBarry Smith MPI_Comm comm; 1000ace3abfcSBarry Smith PetscBool debug = PETSC_FALSE; 100189d82c54SBarry Smith 100289d82c54SBarry Smith PetscFunctionBegin; 1003ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)mapping,&comm);CHKERRQ(ierr); 100424cf384cSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 100524cf384cSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 100624cf384cSBarry Smith if (size == 1) { 100724cf384cSBarry Smith *nproc = 0; 10080298fd71SBarry Smith *procs = NULL; 100995dccacaSBarry Smith ierr = PetscNew(numprocs);CHKERRQ(ierr); 10101e2105dcSBarry Smith (*numprocs)[0] = 0; 101195dccacaSBarry Smith ierr = PetscNew(indices);CHKERRQ(ierr); 10120298fd71SBarry Smith (*indices)[0] = NULL; 1013268a049cSStefano Zampini /* save info for reuse */ 1014268a049cSStefano Zampini mapping->info_nproc = *nproc; 1015268a049cSStefano Zampini mapping->info_procs = *procs; 1016268a049cSStefano Zampini mapping->info_numprocs = *numprocs; 1017268a049cSStefano Zampini mapping->info_indices = *indices; 1018268a049cSStefano Zampini mapping->info_cached = PETSC_TRUE; 101924cf384cSBarry Smith PetscFunctionReturn(0); 102024cf384cSBarry Smith } 102124cf384cSBarry Smith 1022c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)mapping)->options,NULL,"-islocaltoglobalmappinggetinfo_debug",&debug,NULL);CHKERRQ(ierr); 102307b52d57SBarry Smith 10243677ff5aSBarry Smith /* 10256a818285SBarry Smith Notes on ISLocalToGlobalMappingGetBlockInfo 10263677ff5aSBarry Smith 10273677ff5aSBarry Smith globally owned node - the nodes that have been assigned to this processor in global 10283677ff5aSBarry Smith numbering, just for this routine. 10293677ff5aSBarry Smith 10303677ff5aSBarry Smith nontrivial globally owned node - node assigned to this processor that is on a subdomain 10313677ff5aSBarry Smith boundary (i.e. is has more than one local owner) 10323677ff5aSBarry Smith 10333677ff5aSBarry Smith locally owned node - node that exists on this processors subdomain 10343677ff5aSBarry Smith 10353677ff5aSBarry Smith nontrivial locally owned node - node that is not in the interior (i.e. has more than one 10363677ff5aSBarry Smith local subdomain 10373677ff5aSBarry Smith */ 103824cf384cSBarry Smith ierr = PetscObjectGetNewTag((PetscObject)mapping,&tag1);CHKERRQ(ierr); 103924cf384cSBarry Smith ierr = PetscObjectGetNewTag((PetscObject)mapping,&tag2);CHKERRQ(ierr); 104024cf384cSBarry Smith ierr = PetscObjectGetNewTag((PetscObject)mapping,&tag3);CHKERRQ(ierr); 104189d82c54SBarry Smith 104289d82c54SBarry Smith for (i=0; i<n; i++) { 104389d82c54SBarry Smith if (lindices[i] > max) max = lindices[i]; 104489d82c54SBarry Smith } 1045b2566f29SBarry Smith ierr = MPIU_Allreduce(&max,&Ng,1,MPIU_INT,MPI_MAX,comm);CHKERRQ(ierr); 104678058e43SBarry Smith Ng++; 104789d82c54SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 104889d82c54SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1049bc8ff85bSBarry Smith scale = Ng/size + 1; 1050a2e34c3dSBarry Smith ng = scale; if (rank == size-1) ng = Ng - scale*(size-1); ng = PetscMax(1,ng); 1051caba0dd0SBarry Smith rstart = scale*rank; 105289d82c54SBarry Smith 105389d82c54SBarry Smith /* determine ownership ranges of global indices */ 1054785e854fSJed Brown ierr = PetscMalloc1(2*size,&nprocs);CHKERRQ(ierr); 105532dcc486SBarry Smith ierr = PetscMemzero(nprocs,2*size*sizeof(PetscInt));CHKERRQ(ierr); 105689d82c54SBarry Smith 105789d82c54SBarry Smith /* determine owners of each local node */ 1058785e854fSJed Brown ierr = PetscMalloc1(n,&owner);CHKERRQ(ierr); 105989d82c54SBarry Smith for (i=0; i<n; i++) { 10603677ff5aSBarry Smith proc = lindices[i]/scale; /* processor that globally owns this index */ 106127c402fcSBarry Smith nprocs[2*proc+1] = 1; /* processor globally owns at least one of ours */ 10623677ff5aSBarry Smith owner[i] = proc; 106327c402fcSBarry Smith nprocs[2*proc]++; /* count of how many that processor globally owns of ours */ 106489d82c54SBarry Smith } 106527c402fcSBarry Smith nsends = 0; for (i=0; i<size; i++) nsends += nprocs[2*i+1]; 10667904a332SBarry Smith ierr = PetscInfo1(mapping,"Number of global owners for my local data %D\n",nsends);CHKERRQ(ierr); 106789d82c54SBarry Smith 106889d82c54SBarry Smith /* inform other processors of number of messages and max length*/ 106927c402fcSBarry Smith ierr = PetscMaxSum(comm,nprocs,&nmax,&nrecvs);CHKERRQ(ierr); 10707904a332SBarry Smith ierr = PetscInfo1(mapping,"Number of local owners for my global data %D\n",nrecvs);CHKERRQ(ierr); 107189d82c54SBarry Smith 107289d82c54SBarry Smith /* post receives for owned rows */ 1073785e854fSJed Brown ierr = PetscMalloc1((2*nrecvs+1)*(nmax+1),&recvs);CHKERRQ(ierr); 1074854ce69bSBarry Smith ierr = PetscMalloc1(nrecvs+1,&recv_waits);CHKERRQ(ierr); 107589d82c54SBarry Smith for (i=0; i<nrecvs; i++) { 107632dcc486SBarry Smith ierr = MPI_Irecv(recvs+2*nmax*i,2*nmax,MPIU_INT,MPI_ANY_SOURCE,tag1,comm,recv_waits+i);CHKERRQ(ierr); 107789d82c54SBarry Smith } 107889d82c54SBarry Smith 107989d82c54SBarry Smith /* pack messages containing lists of local nodes to owners */ 1080854ce69bSBarry Smith ierr = PetscMalloc1(2*n+1,&sends);CHKERRQ(ierr); 1081854ce69bSBarry Smith ierr = PetscMalloc1(size+1,&starts);CHKERRQ(ierr); 108289d82c54SBarry Smith starts[0] = 0; 1083f6e5521dSKarl Rupp for (i=1; i<size; i++) starts[i] = starts[i-1] + 2*nprocs[2*i-2]; 108489d82c54SBarry Smith for (i=0; i<n; i++) { 108589d82c54SBarry Smith sends[starts[owner[i]]++] = lindices[i]; 108630dcb7c9SBarry Smith sends[starts[owner[i]]++] = i; 108789d82c54SBarry Smith } 108889d82c54SBarry Smith ierr = PetscFree(owner);CHKERRQ(ierr); 108989d82c54SBarry Smith starts[0] = 0; 1090f6e5521dSKarl Rupp for (i=1; i<size; i++) starts[i] = starts[i-1] + 2*nprocs[2*i-2]; 109189d82c54SBarry Smith 109289d82c54SBarry Smith /* send the messages */ 1093854ce69bSBarry Smith ierr = PetscMalloc1(nsends+1,&send_waits);CHKERRQ(ierr); 1094854ce69bSBarry Smith ierr = PetscMalloc1(nsends+1,&dest);CHKERRQ(ierr); 109589d82c54SBarry Smith cnt = 0; 109689d82c54SBarry Smith for (i=0; i<size; i++) { 109727c402fcSBarry Smith if (nprocs[2*i]) { 109832dcc486SBarry Smith ierr = MPI_Isend(sends+starts[i],2*nprocs[2*i],MPIU_INT,i,tag1,comm,send_waits+cnt);CHKERRQ(ierr); 109930dcb7c9SBarry Smith dest[cnt] = i; 110089d82c54SBarry Smith cnt++; 110189d82c54SBarry Smith } 110289d82c54SBarry Smith } 110389d82c54SBarry Smith ierr = PetscFree(starts);CHKERRQ(ierr); 110489d82c54SBarry Smith 110589d82c54SBarry Smith /* wait on receives */ 1106854ce69bSBarry Smith ierr = PetscMalloc1(nrecvs+1,&source);CHKERRQ(ierr); 1107854ce69bSBarry Smith ierr = PetscMalloc1(nrecvs+1,&len);CHKERRQ(ierr); 110889d82c54SBarry Smith cnt = nrecvs; 1109854ce69bSBarry Smith ierr = PetscMalloc1(ng+1,&nownedsenders);CHKERRQ(ierr); 111032dcc486SBarry Smith ierr = PetscMemzero(nownedsenders,ng*sizeof(PetscInt));CHKERRQ(ierr); 111189d82c54SBarry Smith while (cnt) { 111289d82c54SBarry Smith ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr); 111389d82c54SBarry Smith /* unpack receives into our local space */ 111432dcc486SBarry Smith ierr = MPI_Get_count(&recv_status,MPIU_INT,&len[imdex]);CHKERRQ(ierr); 111589d82c54SBarry Smith source[imdex] = recv_status.MPI_SOURCE; 111630dcb7c9SBarry Smith len[imdex] = len[imdex]/2; 1117caba0dd0SBarry Smith /* count how many local owners for each of my global owned indices */ 111830dcb7c9SBarry Smith for (i=0; i<len[imdex]; i++) nownedsenders[recvs[2*imdex*nmax+2*i]-rstart]++; 111989d82c54SBarry Smith cnt--; 112089d82c54SBarry Smith } 112189d82c54SBarry Smith ierr = PetscFree(recv_waits);CHKERRQ(ierr); 112289d82c54SBarry Smith 112330dcb7c9SBarry Smith /* count how many globally owned indices are on an edge multiplied by how many processors own them. */ 1124bc8ff85bSBarry Smith nowned = 0; 1125bc8ff85bSBarry Smith nownedm = 0; 1126bc8ff85bSBarry Smith for (i=0; i<ng; i++) { 1127bc8ff85bSBarry Smith if (nownedsenders[i] > 1) {nownedm += nownedsenders[i]; nowned++;} 1128bc8ff85bSBarry Smith } 1129bc8ff85bSBarry Smith 1130bc8ff85bSBarry Smith /* create single array to contain rank of all local owners of each globally owned index */ 1131854ce69bSBarry Smith ierr = PetscMalloc1(nownedm+1,&ownedsenders);CHKERRQ(ierr); 1132854ce69bSBarry Smith ierr = PetscMalloc1(ng+1,&starts);CHKERRQ(ierr); 1133bc8ff85bSBarry Smith starts[0] = 0; 1134bc8ff85bSBarry Smith for (i=1; i<ng; i++) { 1135bc8ff85bSBarry Smith if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1]; 1136bc8ff85bSBarry Smith else starts[i] = starts[i-1]; 1137bc8ff85bSBarry Smith } 1138bc8ff85bSBarry Smith 113930dcb7c9SBarry Smith /* for each nontrival globally owned node list all arriving processors */ 1140bc8ff85bSBarry Smith for (i=0; i<nrecvs; i++) { 1141bc8ff85bSBarry Smith for (j=0; j<len[i]; j++) { 114230dcb7c9SBarry Smith node = recvs[2*i*nmax+2*j]-rstart; 1143f6e5521dSKarl Rupp if (nownedsenders[node] > 1) ownedsenders[starts[node]++] = source[i]; 1144bc8ff85bSBarry Smith } 1145bc8ff85bSBarry Smith } 1146bc8ff85bSBarry Smith 114707b52d57SBarry Smith if (debug) { /* ----------------------------------- */ 114830dcb7c9SBarry Smith starts[0] = 0; 114930dcb7c9SBarry Smith for (i=1; i<ng; i++) { 115030dcb7c9SBarry Smith if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1]; 115130dcb7c9SBarry Smith else starts[i] = starts[i-1]; 115230dcb7c9SBarry Smith } 115330dcb7c9SBarry Smith for (i=0; i<ng; i++) { 115430dcb7c9SBarry Smith if (nownedsenders[i] > 1) { 11557904a332SBarry Smith ierr = PetscSynchronizedPrintf(comm,"[%d] global node %D local owner processors: ",rank,i+rstart);CHKERRQ(ierr); 115630dcb7c9SBarry Smith for (j=0; j<nownedsenders[i]; j++) { 11577904a332SBarry Smith ierr = PetscSynchronizedPrintf(comm,"%D ",ownedsenders[starts[i]+j]);CHKERRQ(ierr); 115830dcb7c9SBarry Smith } 115930dcb7c9SBarry Smith ierr = PetscSynchronizedPrintf(comm,"\n");CHKERRQ(ierr); 116030dcb7c9SBarry Smith } 116130dcb7c9SBarry Smith } 11620ec8b6e3SBarry Smith ierr = PetscSynchronizedFlush(comm,PETSC_STDOUT);CHKERRQ(ierr); 116307b52d57SBarry Smith } /* ----------------------------------- */ 116430dcb7c9SBarry Smith 11653677ff5aSBarry Smith /* wait on original sends */ 11663a96401aSBarry Smith if (nsends) { 1167785e854fSJed Brown ierr = PetscMalloc1(nsends,&send_status);CHKERRQ(ierr); 11683a96401aSBarry Smith ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr); 11693a96401aSBarry Smith ierr = PetscFree(send_status);CHKERRQ(ierr); 11703a96401aSBarry Smith } 117189d82c54SBarry Smith ierr = PetscFree(send_waits);CHKERRQ(ierr); 11723a96401aSBarry Smith ierr = PetscFree(sends);CHKERRQ(ierr); 11733677ff5aSBarry Smith ierr = PetscFree(nprocs);CHKERRQ(ierr); 11743677ff5aSBarry Smith 11753677ff5aSBarry Smith /* pack messages to send back to local owners */ 117630dcb7c9SBarry Smith starts[0] = 0; 117730dcb7c9SBarry Smith for (i=1; i<ng; i++) { 117830dcb7c9SBarry Smith if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1]; 117930dcb7c9SBarry Smith else starts[i] = starts[i-1]; 118030dcb7c9SBarry Smith } 118130dcb7c9SBarry Smith nsends2 = nrecvs; 1182854ce69bSBarry Smith ierr = PetscMalloc1(nsends2+1,&nprocs);CHKERRQ(ierr); /* length of each message */ 118330dcb7c9SBarry Smith for (i=0; i<nrecvs; i++) { 118430dcb7c9SBarry Smith nprocs[i] = 1; 118530dcb7c9SBarry Smith for (j=0; j<len[i]; j++) { 118630dcb7c9SBarry Smith node = recvs[2*i*nmax+2*j]-rstart; 1187f6e5521dSKarl Rupp if (nownedsenders[node] > 1) nprocs[i] += 2 + nownedsenders[node]; 118830dcb7c9SBarry Smith } 118930dcb7c9SBarry Smith } 1190f6e5521dSKarl Rupp nt = 0; 1191f6e5521dSKarl Rupp for (i=0; i<nsends2; i++) nt += nprocs[i]; 1192f6e5521dSKarl Rupp 1193854ce69bSBarry Smith ierr = PetscMalloc1(nt+1,&sends2);CHKERRQ(ierr); 1194854ce69bSBarry Smith ierr = PetscMalloc1(nsends2+1,&starts2);CHKERRQ(ierr); 1195f6e5521dSKarl Rupp 1196f6e5521dSKarl Rupp starts2[0] = 0; 1197f6e5521dSKarl Rupp for (i=1; i<nsends2; i++) starts2[i] = starts2[i-1] + nprocs[i-1]; 119830dcb7c9SBarry Smith /* 119930dcb7c9SBarry Smith Each message is 1 + nprocs[i] long, and consists of 120030dcb7c9SBarry Smith (0) the number of nodes being sent back 120130dcb7c9SBarry Smith (1) the local node number, 120230dcb7c9SBarry Smith (2) the number of processors sharing it, 120330dcb7c9SBarry Smith (3) the processors sharing it 120430dcb7c9SBarry Smith */ 120530dcb7c9SBarry Smith for (i=0; i<nsends2; i++) { 120630dcb7c9SBarry Smith cnt = 1; 120730dcb7c9SBarry Smith sends2[starts2[i]] = 0; 120830dcb7c9SBarry Smith for (j=0; j<len[i]; j++) { 120930dcb7c9SBarry Smith node = recvs[2*i*nmax+2*j]-rstart; 121030dcb7c9SBarry Smith if (nownedsenders[node] > 1) { 121130dcb7c9SBarry Smith sends2[starts2[i]]++; 121230dcb7c9SBarry Smith sends2[starts2[i]+cnt++] = recvs[2*i*nmax+2*j+1]; 121330dcb7c9SBarry Smith sends2[starts2[i]+cnt++] = nownedsenders[node]; 121432dcc486SBarry Smith ierr = PetscMemcpy(&sends2[starts2[i]+cnt],&ownedsenders[starts[node]],nownedsenders[node]*sizeof(PetscInt));CHKERRQ(ierr); 121530dcb7c9SBarry Smith cnt += nownedsenders[node]; 121630dcb7c9SBarry Smith } 121730dcb7c9SBarry Smith } 121830dcb7c9SBarry Smith } 121930dcb7c9SBarry Smith 122030dcb7c9SBarry Smith /* receive the message lengths */ 122130dcb7c9SBarry Smith nrecvs2 = nsends; 1222854ce69bSBarry Smith ierr = PetscMalloc1(nrecvs2+1,&lens2);CHKERRQ(ierr); 1223854ce69bSBarry Smith ierr = PetscMalloc1(nrecvs2+1,&starts3);CHKERRQ(ierr); 1224854ce69bSBarry Smith ierr = PetscMalloc1(nrecvs2+1,&recv_waits);CHKERRQ(ierr); 122530dcb7c9SBarry Smith for (i=0; i<nrecvs2; i++) { 1226d44834fbSBarry Smith ierr = MPI_Irecv(&lens2[i],1,MPIU_INT,dest[i],tag2,comm,recv_waits+i);CHKERRQ(ierr); 122730dcb7c9SBarry Smith } 1228d44834fbSBarry Smith 12298a8e0b3aSBarry Smith /* send the message lengths */ 12308a8e0b3aSBarry Smith for (i=0; i<nsends2; i++) { 12318a8e0b3aSBarry Smith ierr = MPI_Send(&nprocs[i],1,MPIU_INT,source[i],tag2,comm);CHKERRQ(ierr); 12328a8e0b3aSBarry Smith } 12338a8e0b3aSBarry Smith 1234d44834fbSBarry Smith /* wait on receives of lens */ 12350c468ba9SBarry Smith if (nrecvs2) { 1236785e854fSJed Brown ierr = PetscMalloc1(nrecvs2,&recv_statuses);CHKERRQ(ierr); 1237d44834fbSBarry Smith ierr = MPI_Waitall(nrecvs2,recv_waits,recv_statuses);CHKERRQ(ierr); 1238d44834fbSBarry Smith ierr = PetscFree(recv_statuses);CHKERRQ(ierr); 12390c468ba9SBarry Smith } 1240a2ea699eSBarry Smith ierr = PetscFree(recv_waits);CHKERRQ(ierr); 1241d44834fbSBarry Smith 124230dcb7c9SBarry Smith starts3[0] = 0; 1243d44834fbSBarry Smith nt = 0; 124430dcb7c9SBarry Smith for (i=0; i<nrecvs2-1; i++) { 124530dcb7c9SBarry Smith starts3[i+1] = starts3[i] + lens2[i]; 1246d44834fbSBarry Smith nt += lens2[i]; 124730dcb7c9SBarry Smith } 124876466f69SStefano Zampini if (nrecvs2) nt += lens2[nrecvs2-1]; 1249d44834fbSBarry Smith 1250854ce69bSBarry Smith ierr = PetscMalloc1(nt+1,&recvs2);CHKERRQ(ierr); 1251854ce69bSBarry Smith ierr = PetscMalloc1(nrecvs2+1,&recv_waits);CHKERRQ(ierr); 125252b72c4aSBarry Smith for (i=0; i<nrecvs2; i++) { 125332dcc486SBarry Smith ierr = MPI_Irecv(recvs2+starts3[i],lens2[i],MPIU_INT,dest[i],tag3,comm,recv_waits+i);CHKERRQ(ierr); 125430dcb7c9SBarry Smith } 125530dcb7c9SBarry Smith 125630dcb7c9SBarry Smith /* send the messages */ 1257854ce69bSBarry Smith ierr = PetscMalloc1(nsends2+1,&send_waits);CHKERRQ(ierr); 125830dcb7c9SBarry Smith for (i=0; i<nsends2; i++) { 125932dcc486SBarry Smith ierr = MPI_Isend(sends2+starts2[i],nprocs[i],MPIU_INT,source[i],tag3,comm,send_waits+i);CHKERRQ(ierr); 126030dcb7c9SBarry Smith } 126130dcb7c9SBarry Smith 126230dcb7c9SBarry Smith /* wait on receives */ 12630c468ba9SBarry Smith if (nrecvs2) { 1264785e854fSJed Brown ierr = PetscMalloc1(nrecvs2,&recv_statuses);CHKERRQ(ierr); 126530dcb7c9SBarry Smith ierr = MPI_Waitall(nrecvs2,recv_waits,recv_statuses);CHKERRQ(ierr); 126630dcb7c9SBarry Smith ierr = PetscFree(recv_statuses);CHKERRQ(ierr); 12670c468ba9SBarry Smith } 126830dcb7c9SBarry Smith ierr = PetscFree(recv_waits);CHKERRQ(ierr); 126930dcb7c9SBarry Smith ierr = PetscFree(nprocs);CHKERRQ(ierr); 127030dcb7c9SBarry Smith 127107b52d57SBarry Smith if (debug) { /* ----------------------------------- */ 127230dcb7c9SBarry Smith cnt = 0; 127330dcb7c9SBarry Smith for (i=0; i<nrecvs2; i++) { 127430dcb7c9SBarry Smith nt = recvs2[cnt++]; 127530dcb7c9SBarry Smith for (j=0; j<nt; j++) { 12767904a332SBarry Smith ierr = PetscSynchronizedPrintf(comm,"[%d] local node %D number of subdomains %D: ",rank,recvs2[cnt],recvs2[cnt+1]);CHKERRQ(ierr); 127730dcb7c9SBarry Smith for (k=0; k<recvs2[cnt+1]; k++) { 12787904a332SBarry Smith ierr = PetscSynchronizedPrintf(comm,"%D ",recvs2[cnt+2+k]);CHKERRQ(ierr); 127930dcb7c9SBarry Smith } 128030dcb7c9SBarry Smith cnt += 2 + recvs2[cnt+1]; 128130dcb7c9SBarry Smith ierr = PetscSynchronizedPrintf(comm,"\n");CHKERRQ(ierr); 128230dcb7c9SBarry Smith } 128330dcb7c9SBarry Smith } 12840ec8b6e3SBarry Smith ierr = PetscSynchronizedFlush(comm,PETSC_STDOUT);CHKERRQ(ierr); 128507b52d57SBarry Smith } /* ----------------------------------- */ 128630dcb7c9SBarry Smith 128730dcb7c9SBarry Smith /* count number subdomains for each local node */ 1288785e854fSJed Brown ierr = PetscMalloc1(size,&nprocs);CHKERRQ(ierr); 128932dcc486SBarry Smith ierr = PetscMemzero(nprocs,size*sizeof(PetscInt));CHKERRQ(ierr); 129030dcb7c9SBarry Smith cnt = 0; 129130dcb7c9SBarry Smith for (i=0; i<nrecvs2; i++) { 129230dcb7c9SBarry Smith nt = recvs2[cnt++]; 129330dcb7c9SBarry Smith for (j=0; j<nt; j++) { 1294f6e5521dSKarl Rupp for (k=0; k<recvs2[cnt+1]; k++) nprocs[recvs2[cnt+2+k]]++; 129530dcb7c9SBarry Smith cnt += 2 + recvs2[cnt+1]; 129630dcb7c9SBarry Smith } 129730dcb7c9SBarry Smith } 129830dcb7c9SBarry Smith nt = 0; for (i=0; i<size; i++) nt += (nprocs[i] > 0); 129930dcb7c9SBarry Smith *nproc = nt; 1300854ce69bSBarry Smith ierr = PetscMalloc1(nt+1,procs);CHKERRQ(ierr); 1301854ce69bSBarry Smith ierr = PetscMalloc1(nt+1,numprocs);CHKERRQ(ierr); 1302854ce69bSBarry Smith ierr = PetscMalloc1(nt+1,indices);CHKERRQ(ierr); 13030298fd71SBarry Smith for (i=0;i<nt+1;i++) (*indices)[i]=NULL; 1304785e854fSJed Brown ierr = PetscMalloc1(size,&bprocs);CHKERRQ(ierr); 130530dcb7c9SBarry Smith cnt = 0; 130630dcb7c9SBarry Smith for (i=0; i<size; i++) { 130730dcb7c9SBarry Smith if (nprocs[i] > 0) { 130830dcb7c9SBarry Smith bprocs[i] = cnt; 130930dcb7c9SBarry Smith (*procs)[cnt] = i; 131030dcb7c9SBarry Smith (*numprocs)[cnt] = nprocs[i]; 1311785e854fSJed Brown ierr = PetscMalloc1(nprocs[i],&(*indices)[cnt]);CHKERRQ(ierr); 131230dcb7c9SBarry Smith cnt++; 131330dcb7c9SBarry Smith } 131430dcb7c9SBarry Smith } 131530dcb7c9SBarry Smith 131630dcb7c9SBarry Smith /* make the list of subdomains for each nontrivial local node */ 131732dcc486SBarry Smith ierr = PetscMemzero(*numprocs,nt*sizeof(PetscInt));CHKERRQ(ierr); 131830dcb7c9SBarry Smith cnt = 0; 131930dcb7c9SBarry Smith for (i=0; i<nrecvs2; i++) { 132030dcb7c9SBarry Smith nt = recvs2[cnt++]; 132130dcb7c9SBarry Smith for (j=0; j<nt; j++) { 1322f6e5521dSKarl Rupp for (k=0; k<recvs2[cnt+1]; k++) (*indices)[bprocs[recvs2[cnt+2+k]]][(*numprocs)[bprocs[recvs2[cnt+2+k]]]++] = recvs2[cnt]; 132330dcb7c9SBarry Smith cnt += 2 + recvs2[cnt+1]; 132430dcb7c9SBarry Smith } 132530dcb7c9SBarry Smith } 132630dcb7c9SBarry Smith ierr = PetscFree(bprocs);CHKERRQ(ierr); 132707b52d57SBarry Smith ierr = PetscFree(recvs2);CHKERRQ(ierr); 132830dcb7c9SBarry Smith 132907b52d57SBarry Smith /* sort the node indexing by their global numbers */ 133007b52d57SBarry Smith nt = *nproc; 133107b52d57SBarry Smith for (i=0; i<nt; i++) { 1332854ce69bSBarry Smith ierr = PetscMalloc1((*numprocs)[i],&tmp);CHKERRQ(ierr); 1333f6e5521dSKarl Rupp for (j=0; j<(*numprocs)[i]; j++) tmp[j] = lindices[(*indices)[i][j]]; 133407b52d57SBarry Smith ierr = PetscSortIntWithArray((*numprocs)[i],tmp,(*indices)[i]);CHKERRQ(ierr); 133507b52d57SBarry Smith ierr = PetscFree(tmp);CHKERRQ(ierr); 133607b52d57SBarry Smith } 133707b52d57SBarry Smith 133807b52d57SBarry Smith if (debug) { /* ----------------------------------- */ 133930dcb7c9SBarry Smith nt = *nproc; 134030dcb7c9SBarry Smith for (i=0; i<nt; i++) { 13417904a332SBarry Smith ierr = PetscSynchronizedPrintf(comm,"[%d] subdomain %D number of indices %D: ",rank,(*procs)[i],(*numprocs)[i]);CHKERRQ(ierr); 134230dcb7c9SBarry Smith for (j=0; j<(*numprocs)[i]; j++) { 13437904a332SBarry Smith ierr = PetscSynchronizedPrintf(comm,"%D ",(*indices)[i][j]);CHKERRQ(ierr); 134430dcb7c9SBarry Smith } 134530dcb7c9SBarry Smith ierr = PetscSynchronizedPrintf(comm,"\n");CHKERRQ(ierr); 134630dcb7c9SBarry Smith } 13470ec8b6e3SBarry Smith ierr = PetscSynchronizedFlush(comm,PETSC_STDOUT);CHKERRQ(ierr); 134807b52d57SBarry Smith } /* ----------------------------------- */ 134930dcb7c9SBarry Smith 135030dcb7c9SBarry Smith /* wait on sends */ 135130dcb7c9SBarry Smith if (nsends2) { 1352785e854fSJed Brown ierr = PetscMalloc1(nsends2,&send_status);CHKERRQ(ierr); 135330dcb7c9SBarry Smith ierr = MPI_Waitall(nsends2,send_waits,send_status);CHKERRQ(ierr); 135430dcb7c9SBarry Smith ierr = PetscFree(send_status);CHKERRQ(ierr); 135530dcb7c9SBarry Smith } 135630dcb7c9SBarry Smith 135730dcb7c9SBarry Smith ierr = PetscFree(starts3);CHKERRQ(ierr); 135830dcb7c9SBarry Smith ierr = PetscFree(dest);CHKERRQ(ierr); 135930dcb7c9SBarry Smith ierr = PetscFree(send_waits);CHKERRQ(ierr); 13603677ff5aSBarry Smith 1361bc8ff85bSBarry Smith ierr = PetscFree(nownedsenders);CHKERRQ(ierr); 1362bc8ff85bSBarry Smith ierr = PetscFree(ownedsenders);CHKERRQ(ierr); 1363bc8ff85bSBarry Smith ierr = PetscFree(starts);CHKERRQ(ierr); 136430dcb7c9SBarry Smith ierr = PetscFree(starts2);CHKERRQ(ierr); 136530dcb7c9SBarry Smith ierr = PetscFree(lens2);CHKERRQ(ierr); 136689d82c54SBarry Smith 136789d82c54SBarry Smith ierr = PetscFree(source);CHKERRQ(ierr); 136897f1f81fSBarry Smith ierr = PetscFree(len);CHKERRQ(ierr); 136989d82c54SBarry Smith ierr = PetscFree(recvs);CHKERRQ(ierr); 13703a96401aSBarry Smith ierr = PetscFree(nprocs);CHKERRQ(ierr); 137130dcb7c9SBarry Smith ierr = PetscFree(sends2);CHKERRQ(ierr); 137224cf384cSBarry Smith 137324cf384cSBarry Smith /* put the information about myself as the first entry in the list */ 137424cf384cSBarry Smith first_procs = (*procs)[0]; 137524cf384cSBarry Smith first_numprocs = (*numprocs)[0]; 137624cf384cSBarry Smith first_indices = (*indices)[0]; 137724cf384cSBarry Smith for (i=0; i<*nproc; i++) { 137824cf384cSBarry Smith if ((*procs)[i] == rank) { 137924cf384cSBarry Smith (*procs)[0] = (*procs)[i]; 138024cf384cSBarry Smith (*numprocs)[0] = (*numprocs)[i]; 138124cf384cSBarry Smith (*indices)[0] = (*indices)[i]; 138224cf384cSBarry Smith (*procs)[i] = first_procs; 138324cf384cSBarry Smith (*numprocs)[i] = first_numprocs; 138424cf384cSBarry Smith (*indices)[i] = first_indices; 138524cf384cSBarry Smith break; 138624cf384cSBarry Smith } 138724cf384cSBarry Smith } 1388268a049cSStefano Zampini 1389268a049cSStefano Zampini /* save info for reuse */ 1390268a049cSStefano Zampini mapping->info_nproc = *nproc; 1391268a049cSStefano Zampini mapping->info_procs = *procs; 1392268a049cSStefano Zampini mapping->info_numprocs = *numprocs; 1393268a049cSStefano Zampini mapping->info_indices = *indices; 1394268a049cSStefano Zampini mapping->info_cached = PETSC_TRUE; 139589d82c54SBarry Smith PetscFunctionReturn(0); 139689d82c54SBarry Smith } 139789d82c54SBarry Smith 13986a818285SBarry Smith /*@C 13996a818285SBarry Smith ISLocalToGlobalMappingRestoreBlockInfo - Frees the memory allocated by ISLocalToGlobalMappingGetBlockInfo() 14006a818285SBarry Smith 14016a818285SBarry Smith Collective on ISLocalToGlobalMapping 14026a818285SBarry Smith 14036a818285SBarry Smith Input Parameters: 14046a818285SBarry Smith . mapping - the mapping from local to global indexing 14056a818285SBarry Smith 14066a818285SBarry Smith Output Parameter: 14076a818285SBarry Smith + nproc - number of processors that are connected to this one 14086a818285SBarry Smith . proc - neighboring processors 14096a818285SBarry Smith . numproc - number of indices for each processor 14106a818285SBarry Smith - indices - indices of local nodes shared with neighbor (sorted by global numbering) 14116a818285SBarry Smith 14126a818285SBarry Smith Level: advanced 14136a818285SBarry Smith 14146a818285SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(), 14156a818285SBarry Smith ISLocalToGlobalMappingGetInfo() 14166a818285SBarry Smith @*/ 14176a818285SBarry Smith PetscErrorCode ISLocalToGlobalMappingRestoreBlockInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[]) 14186a818285SBarry Smith { 14196a818285SBarry Smith PetscErrorCode ierr; 14206a818285SBarry Smith 14216a818285SBarry Smith PetscFunctionBegin; 1422cbc1caf0SMatthew G. Knepley PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1); 1423268a049cSStefano Zampini if (mapping->info_free) { 14246a818285SBarry Smith ierr = PetscFree(*numprocs);CHKERRQ(ierr); 14256a818285SBarry Smith if (*indices) { 1426268a049cSStefano Zampini PetscInt i; 1427268a049cSStefano Zampini 14286a818285SBarry Smith ierr = PetscFree((*indices)[0]);CHKERRQ(ierr); 14296a818285SBarry Smith for (i=1; i<*nproc; i++) { 14306a818285SBarry Smith ierr = PetscFree((*indices)[i]);CHKERRQ(ierr); 14316a818285SBarry Smith } 14326a818285SBarry Smith ierr = PetscFree(*indices);CHKERRQ(ierr); 14336a818285SBarry Smith } 1434268a049cSStefano Zampini } 1435268a049cSStefano Zampini *nproc = 0; 1436268a049cSStefano Zampini *procs = NULL; 1437268a049cSStefano Zampini *numprocs = NULL; 1438268a049cSStefano Zampini *indices = NULL; 14396a818285SBarry Smith PetscFunctionReturn(0); 14406a818285SBarry Smith } 14416a818285SBarry Smith 14426a818285SBarry Smith /*@C 14436a818285SBarry Smith ISLocalToGlobalMappingGetInfo - Gets the neighbor information for each processor and 14446a818285SBarry Smith each index shared by more than one processor 14456a818285SBarry Smith 14466a818285SBarry Smith Collective on ISLocalToGlobalMapping 14476a818285SBarry Smith 14486a818285SBarry Smith Input Parameters: 14496a818285SBarry Smith . mapping - the mapping from local to global indexing 14506a818285SBarry Smith 14516a818285SBarry Smith Output Parameter: 14526a818285SBarry Smith + nproc - number of processors that are connected to this one 14536a818285SBarry Smith . proc - neighboring processors 14546a818285SBarry Smith . numproc - number of indices for each subdomain (processor) 14556a818285SBarry Smith - indices - indices of nodes (in local numbering) shared with neighbors (sorted by global numbering) 14566a818285SBarry Smith 14576a818285SBarry Smith Level: advanced 14586a818285SBarry Smith 14596a818285SBarry Smith Concepts: mapping^local to global 14606a818285SBarry Smith 14611bd0b88eSStefano Zampini Notes: The user needs to call ISLocalToGlobalMappingRestoreInfo when the data is no longer needed. 14621bd0b88eSStefano Zampini 14636a818285SBarry Smith Fortran Usage: 14646a818285SBarry Smith $ ISLocalToGlobalMpngGetInfoSize(ISLocalToGlobalMapping,PetscInt nproc,PetscInt numprocmax,ierr) followed by 14656a818285SBarry Smith $ ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt nproc, PetscInt procs[nproc],PetscInt numprocs[nproc], 14666a818285SBarry Smith PetscInt indices[nproc][numprocmax],ierr) 14676a818285SBarry Smith There is no ISLocalToGlobalMappingRestoreInfo() in Fortran. You must make sure that procs[], numprocs[] and 14686a818285SBarry Smith indices[][] are large enough arrays, either by allocating them dynamically or defining static ones large enough. 14696a818285SBarry Smith 14706a818285SBarry Smith 14716a818285SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(), 14726a818285SBarry Smith ISLocalToGlobalMappingRestoreInfo() 14736a818285SBarry Smith @*/ 14746a818285SBarry Smith PetscErrorCode ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[]) 14756a818285SBarry Smith { 14766a818285SBarry Smith PetscErrorCode ierr; 1477268a049cSStefano Zampini PetscInt **bindices = NULL,*bnumprocs = NULL,bs = mapping->bs,i,j,k; 14786a818285SBarry Smith 14796a818285SBarry Smith PetscFunctionBegin; 14806a818285SBarry Smith PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1); 1481268a049cSStefano Zampini ierr = ISLocalToGlobalMappingGetBlockInfo(mapping,nproc,procs,&bnumprocs,&bindices);CHKERRQ(ierr); 1482268a049cSStefano Zampini if (bs > 1) { /* we need to expand the cached info */ 1483732f65e3SBarry Smith ierr = PetscCalloc1(*nproc,&*indices);CHKERRQ(ierr); 1484268a049cSStefano Zampini ierr = PetscCalloc1(*nproc,&*numprocs);CHKERRQ(ierr); 14856a818285SBarry Smith for (i=0; i<*nproc; i++) { 1486268a049cSStefano Zampini ierr = PetscMalloc1(bs*bnumprocs[i],&(*indices)[i]);CHKERRQ(ierr); 1487268a049cSStefano Zampini for (j=0; j<bnumprocs[i]; j++) { 14886a818285SBarry Smith for (k=0; k<bs; k++) { 14896a818285SBarry Smith (*indices)[i][j*bs+k] = bs*bindices[i][j] + k; 14906a818285SBarry Smith } 14916a818285SBarry Smith } 1492268a049cSStefano Zampini (*numprocs)[i] = bnumprocs[i]*bs; 14936a818285SBarry Smith } 1494268a049cSStefano Zampini mapping->info_free = PETSC_TRUE; 1495268a049cSStefano Zampini } else { 1496268a049cSStefano Zampini *numprocs = bnumprocs; 1497268a049cSStefano Zampini *indices = bindices; 14986a818285SBarry Smith } 14996a818285SBarry Smith PetscFunctionReturn(0); 15006a818285SBarry Smith } 15016a818285SBarry Smith 150207b52d57SBarry Smith /*@C 150307b52d57SBarry Smith ISLocalToGlobalMappingRestoreInfo - Frees the memory allocated by ISLocalToGlobalMappingGetInfo() 150489d82c54SBarry Smith 150507b52d57SBarry Smith Collective on ISLocalToGlobalMapping 150607b52d57SBarry Smith 150707b52d57SBarry Smith Input Parameters: 150807b52d57SBarry Smith . mapping - the mapping from local to global indexing 150907b52d57SBarry Smith 151007b52d57SBarry Smith Output Parameter: 151107b52d57SBarry Smith + nproc - number of processors that are connected to this one 151207b52d57SBarry Smith . proc - neighboring processors 151307b52d57SBarry Smith . numproc - number of indices for each processor 151407b52d57SBarry Smith - indices - indices of local nodes shared with neighbor (sorted by global numbering) 151507b52d57SBarry Smith 151607b52d57SBarry Smith Level: advanced 151707b52d57SBarry Smith 151807b52d57SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(), 151907b52d57SBarry Smith ISLocalToGlobalMappingGetInfo() 152007b52d57SBarry Smith @*/ 15217087cfbeSBarry Smith PetscErrorCode ISLocalToGlobalMappingRestoreInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[]) 152207b52d57SBarry Smith { 15236849ba73SBarry Smith PetscErrorCode ierr; 152407b52d57SBarry Smith 152507b52d57SBarry Smith PetscFunctionBegin; 15266a818285SBarry Smith ierr = ISLocalToGlobalMappingRestoreBlockInfo(mapping,nproc,procs,numprocs,indices);CHKERRQ(ierr); 152707b52d57SBarry Smith PetscFunctionReturn(0); 152807b52d57SBarry Smith } 152986994e45SJed Brown 153086994e45SJed Brown /*@C 15311bd0b88eSStefano Zampini ISLocalToGlobalMappingGetNodeInfo - Gets the neighbor information for each node 15321bd0b88eSStefano Zampini 15331bd0b88eSStefano Zampini Collective on ISLocalToGlobalMapping 15341bd0b88eSStefano Zampini 15351bd0b88eSStefano Zampini Input Parameters: 15361bd0b88eSStefano Zampini . mapping - the mapping from local to global indexing 15371bd0b88eSStefano Zampini 15381bd0b88eSStefano Zampini Output Parameter: 15391bd0b88eSStefano Zampini + nnodes - number of local nodes (same ISLocalToGlobalMappingGetSize()) 15401bd0b88eSStefano Zampini . count - number of neighboring processors per node 15411bd0b88eSStefano Zampini - indices - indices of processes sharing the node (sorted) 15421bd0b88eSStefano Zampini 15431bd0b88eSStefano Zampini Level: advanced 15441bd0b88eSStefano Zampini 15451bd0b88eSStefano Zampini Concepts: mapping^local to global 15461bd0b88eSStefano Zampini 15471bd0b88eSStefano Zampini Notes: The user needs to call ISLocalToGlobalMappingRestoreInfo when the data is no longer needed. 15481bd0b88eSStefano Zampini 15491bd0b88eSStefano Zampini .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(), 15501bd0b88eSStefano Zampini ISLocalToGlobalMappingGetInfo(), ISLocalToGlobalMappingRestoreNodeInfo() 15511bd0b88eSStefano Zampini @*/ 15521bd0b88eSStefano Zampini PetscErrorCode ISLocalToGlobalMappingGetNodeInfo(ISLocalToGlobalMapping mapping,PetscInt *nnodes,PetscInt *count[],PetscInt **indices[]) 15531bd0b88eSStefano Zampini { 15541bd0b88eSStefano Zampini PetscInt n; 15551bd0b88eSStefano Zampini PetscErrorCode ierr; 15561bd0b88eSStefano Zampini 15571bd0b88eSStefano Zampini PetscFunctionBegin; 15581bd0b88eSStefano Zampini PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1); 15591bd0b88eSStefano Zampini ierr = ISLocalToGlobalMappingGetSize(mapping,&n);CHKERRQ(ierr); 15601bd0b88eSStefano Zampini if (!mapping->info_nodec) { 15611bd0b88eSStefano Zampini PetscInt i,m,n_neigh,*neigh,*n_shared,**shared; 15621bd0b88eSStefano Zampini 15631bd0b88eSStefano Zampini ierr = PetscCalloc1(n+1,&mapping->info_nodec);CHKERRQ(ierr); /* always allocate to flag setup */ 15641bd0b88eSStefano Zampini ierr = PetscMalloc1(n,&mapping->info_nodei);CHKERRQ(ierr); 15651bd0b88eSStefano Zampini ierr = ISLocalToGlobalMappingGetInfo(mapping,&n_neigh,&neigh,&n_shared,&shared);CHKERRQ(ierr); 15661bd0b88eSStefano Zampini for (i=0,m=0;i<n;i++) { mapping->info_nodec[i] = 1; m++; } 15671bd0b88eSStefano Zampini for (i=1;i<n_neigh;i++) { 15681bd0b88eSStefano Zampini PetscInt j; 15691bd0b88eSStefano Zampini 15701bd0b88eSStefano Zampini m += n_shared[i]; 15711bd0b88eSStefano Zampini for (j=0;j<n_shared[i];j++) mapping->info_nodec[shared[i][j]] += 1; 15721bd0b88eSStefano Zampini } 15731bd0b88eSStefano Zampini if (n) { ierr = PetscMalloc1(m,&mapping->info_nodei[0]);CHKERRQ(ierr); } 15741bd0b88eSStefano Zampini for (i=1;i<n;i++) mapping->info_nodei[i] = mapping->info_nodei[i-1] + mapping->info_nodec[i-1]; 15751bd0b88eSStefano Zampini ierr = PetscMemzero(mapping->info_nodec,n*sizeof(PetscInt));CHKERRQ(ierr); 15761bd0b88eSStefano Zampini for (i=0;i<n;i++) { mapping->info_nodec[i] = 1; mapping->info_nodei[i][0] = neigh[0]; } 15771bd0b88eSStefano Zampini for (i=1;i<n_neigh;i++) { 15781bd0b88eSStefano Zampini PetscInt j; 15791bd0b88eSStefano Zampini 15801bd0b88eSStefano Zampini for (j=0;j<n_shared[i];j++) { 15811bd0b88eSStefano Zampini PetscInt k = shared[i][j]; 15821bd0b88eSStefano Zampini 15831bd0b88eSStefano Zampini mapping->info_nodei[k][mapping->info_nodec[k]] = neigh[i]; 15841bd0b88eSStefano Zampini mapping->info_nodec[k] += 1; 15851bd0b88eSStefano Zampini } 15861bd0b88eSStefano Zampini } 15871bd0b88eSStefano Zampini for (i=0;i<n;i++) { ierr = PetscSortRemoveDupsInt(&mapping->info_nodec[i],mapping->info_nodei[i]);CHKERRQ(ierr); } 15881bd0b88eSStefano Zampini ierr = ISLocalToGlobalMappingRestoreInfo(mapping,&n_neigh,&neigh,&n_shared,&shared);CHKERRQ(ierr); 15891bd0b88eSStefano Zampini } 15901bd0b88eSStefano Zampini if (nnodes) *nnodes = n; 15911bd0b88eSStefano Zampini if (count) *count = mapping->info_nodec; 15921bd0b88eSStefano Zampini if (indices) *indices = mapping->info_nodei; 15931bd0b88eSStefano Zampini PetscFunctionReturn(0); 15941bd0b88eSStefano Zampini } 15951bd0b88eSStefano Zampini 15961bd0b88eSStefano Zampini /*@C 15971bd0b88eSStefano Zampini ISLocalToGlobalMappingRestoreNodeInfo - Frees the memory allocated by ISLocalToGlobalMappingGetNodeInfo() 15981bd0b88eSStefano Zampini 15991bd0b88eSStefano Zampini Collective on ISLocalToGlobalMapping 16001bd0b88eSStefano Zampini 16011bd0b88eSStefano Zampini Input Parameters: 16021bd0b88eSStefano Zampini . mapping - the mapping from local to global indexing 16031bd0b88eSStefano Zampini 16041bd0b88eSStefano Zampini Output Parameter: 16051bd0b88eSStefano Zampini + nnodes - number of local nodes 16061bd0b88eSStefano Zampini . count - number of neighboring processors per node 16071bd0b88eSStefano Zampini - indices - indices of processes sharing the node (sorted) 16081bd0b88eSStefano Zampini 16091bd0b88eSStefano Zampini Level: advanced 16101bd0b88eSStefano Zampini 16111bd0b88eSStefano Zampini .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(), 16121bd0b88eSStefano Zampini ISLocalToGlobalMappingGetInfo() 16131bd0b88eSStefano Zampini @*/ 16141bd0b88eSStefano Zampini PetscErrorCode ISLocalToGlobalMappingRestoreNodeInfo(ISLocalToGlobalMapping mapping,PetscInt *nnodes,PetscInt *count[],PetscInt **indices[]) 16151bd0b88eSStefano Zampini { 16161bd0b88eSStefano Zampini PetscFunctionBegin; 16171bd0b88eSStefano Zampini PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1); 16181bd0b88eSStefano Zampini if (nnodes) *nnodes = 0; 16191bd0b88eSStefano Zampini if (count) *count = NULL; 16201bd0b88eSStefano Zampini if (indices) *indices = NULL; 16211bd0b88eSStefano Zampini PetscFunctionReturn(0); 16221bd0b88eSStefano Zampini } 16231bd0b88eSStefano Zampini 16241bd0b88eSStefano Zampini /*@C 1625107e9a97SBarry Smith ISLocalToGlobalMappingGetIndices - Get global indices for every local point that is mapped 162686994e45SJed Brown 162786994e45SJed Brown Not Collective 162886994e45SJed Brown 162986994e45SJed Brown Input Arguments: 163086994e45SJed Brown . ltog - local to global mapping 163186994e45SJed Brown 163286994e45SJed Brown Output Arguments: 1633565245c5SBarry Smith . array - array of indices, the length of this array may be obtained with ISLocalToGlobalMappingGetSize() 163486994e45SJed Brown 163586994e45SJed Brown Level: advanced 163686994e45SJed Brown 163795452b02SPatrick Sanan Notes: 163895452b02SPatrick Sanan ISLocalToGlobalMappingGetSize() returns the length the this array 1639107e9a97SBarry Smith 1640107e9a97SBarry Smith .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingRestoreIndices(), ISLocalToGlobalMappingGetBlockIndices(), ISLocalToGlobalMappingRestoreBlockIndices() 164186994e45SJed Brown @*/ 16427087cfbeSBarry Smith PetscErrorCode ISLocalToGlobalMappingGetIndices(ISLocalToGlobalMapping ltog,const PetscInt **array) 164386994e45SJed Brown { 164486994e45SJed Brown PetscFunctionBegin; 164586994e45SJed Brown PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1); 164686994e45SJed Brown PetscValidPointer(array,2); 164745b6f7e9SBarry Smith if (ltog->bs == 1) { 164886994e45SJed Brown *array = ltog->indices; 164945b6f7e9SBarry Smith } else { 165045b6f7e9SBarry Smith PetscInt *jj,k,i,j,n = ltog->n, bs = ltog->bs; 165145b6f7e9SBarry Smith const PetscInt *ii; 165245b6f7e9SBarry Smith PetscErrorCode ierr; 165345b6f7e9SBarry Smith 165445b6f7e9SBarry Smith ierr = PetscMalloc1(bs*n,&jj);CHKERRQ(ierr); 165545b6f7e9SBarry Smith *array = jj; 165645b6f7e9SBarry Smith k = 0; 165745b6f7e9SBarry Smith ii = ltog->indices; 165845b6f7e9SBarry Smith for (i=0; i<n; i++) 165945b6f7e9SBarry Smith for (j=0; j<bs; j++) 166045b6f7e9SBarry Smith jj[k++] = bs*ii[i] + j; 166145b6f7e9SBarry Smith } 166286994e45SJed Brown PetscFunctionReturn(0); 166386994e45SJed Brown } 166486994e45SJed Brown 166586994e45SJed Brown /*@C 1666193a2b41SJulian Andrej ISLocalToGlobalMappingRestoreIndices - Restore indices obtained with ISLocalToGlobalMappingGetIndices() 166786994e45SJed Brown 166886994e45SJed Brown Not Collective 166986994e45SJed Brown 167086994e45SJed Brown Input Arguments: 167186994e45SJed Brown + ltog - local to global mapping 167286994e45SJed Brown - array - array of indices 167386994e45SJed Brown 167486994e45SJed Brown Level: advanced 167586994e45SJed Brown 167686994e45SJed Brown .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingGetIndices() 167786994e45SJed Brown @*/ 16787087cfbeSBarry Smith PetscErrorCode ISLocalToGlobalMappingRestoreIndices(ISLocalToGlobalMapping ltog,const PetscInt **array) 167986994e45SJed Brown { 168086994e45SJed Brown PetscFunctionBegin; 168186994e45SJed Brown PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1); 168286994e45SJed Brown PetscValidPointer(array,2); 168345b6f7e9SBarry Smith if (ltog->bs == 1 && *array != ltog->indices) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Trying to return mismatched pointer"); 168445b6f7e9SBarry Smith 168545b6f7e9SBarry Smith if (ltog->bs > 1) { 168645b6f7e9SBarry Smith PetscErrorCode ierr; 168745b6f7e9SBarry Smith ierr = PetscFree(*(void**)array);CHKERRQ(ierr); 168845b6f7e9SBarry Smith } 168945b6f7e9SBarry Smith PetscFunctionReturn(0); 169045b6f7e9SBarry Smith } 169145b6f7e9SBarry Smith 169245b6f7e9SBarry Smith /*@C 169345b6f7e9SBarry Smith ISLocalToGlobalMappingGetBlockIndices - Get global indices for every local block 169445b6f7e9SBarry Smith 169545b6f7e9SBarry Smith Not Collective 169645b6f7e9SBarry Smith 169745b6f7e9SBarry Smith Input Arguments: 169845b6f7e9SBarry Smith . ltog - local to global mapping 169945b6f7e9SBarry Smith 170045b6f7e9SBarry Smith Output Arguments: 170145b6f7e9SBarry Smith . array - array of indices 170245b6f7e9SBarry Smith 170345b6f7e9SBarry Smith Level: advanced 170445b6f7e9SBarry Smith 170545b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingRestoreBlockIndices() 170645b6f7e9SBarry Smith @*/ 170745b6f7e9SBarry Smith PetscErrorCode ISLocalToGlobalMappingGetBlockIndices(ISLocalToGlobalMapping ltog,const PetscInt **array) 170845b6f7e9SBarry Smith { 170945b6f7e9SBarry Smith PetscFunctionBegin; 171045b6f7e9SBarry Smith PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1); 171145b6f7e9SBarry Smith PetscValidPointer(array,2); 171245b6f7e9SBarry Smith *array = ltog->indices; 171345b6f7e9SBarry Smith PetscFunctionReturn(0); 171445b6f7e9SBarry Smith } 171545b6f7e9SBarry Smith 171645b6f7e9SBarry Smith /*@C 171745b6f7e9SBarry Smith ISLocalToGlobalMappingRestoreBlockIndices - Restore indices obtained with ISLocalToGlobalMappingGetBlockIndices() 171845b6f7e9SBarry Smith 171945b6f7e9SBarry Smith Not Collective 172045b6f7e9SBarry Smith 172145b6f7e9SBarry Smith Input Arguments: 172245b6f7e9SBarry Smith + ltog - local to global mapping 172345b6f7e9SBarry Smith - array - array of indices 172445b6f7e9SBarry Smith 172545b6f7e9SBarry Smith Level: advanced 172645b6f7e9SBarry Smith 172745b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingGetIndices() 172845b6f7e9SBarry Smith @*/ 172945b6f7e9SBarry Smith PetscErrorCode ISLocalToGlobalMappingRestoreBlockIndices(ISLocalToGlobalMapping ltog,const PetscInt **array) 173045b6f7e9SBarry Smith { 173145b6f7e9SBarry Smith PetscFunctionBegin; 173245b6f7e9SBarry Smith PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1); 173345b6f7e9SBarry Smith PetscValidPointer(array,2); 173486994e45SJed Brown if (*array != ltog->indices) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Trying to return mismatched pointer"); 17350298fd71SBarry Smith *array = NULL; 173686994e45SJed Brown PetscFunctionReturn(0); 173786994e45SJed Brown } 1738f7efa3c7SJed Brown 1739f7efa3c7SJed Brown /*@C 1740f7efa3c7SJed Brown ISLocalToGlobalMappingConcatenate - Create a new mapping that concatenates a list of mappings 1741f7efa3c7SJed Brown 1742f7efa3c7SJed Brown Not Collective 1743f7efa3c7SJed Brown 1744f7efa3c7SJed Brown Input Arguments: 1745f7efa3c7SJed Brown + comm - communicator for the new mapping, must contain the communicator of every mapping to concatenate 1746f7efa3c7SJed Brown . n - number of mappings to concatenate 1747f7efa3c7SJed Brown - ltogs - local to global mappings 1748f7efa3c7SJed Brown 1749f7efa3c7SJed Brown Output Arguments: 1750f7efa3c7SJed Brown . ltogcat - new mapping 1751f7efa3c7SJed Brown 17529d90f715SBarry Smith Note: this currently always returns a mapping with block size of 1 17539d90f715SBarry Smith 17549d90f715SBarry Smith Developer Note: If all the input mapping have the same block size we could easily handle that as a special case 17559d90f715SBarry Smith 1756f7efa3c7SJed Brown Level: advanced 1757f7efa3c7SJed Brown 1758f7efa3c7SJed Brown .seealso: ISLocalToGlobalMappingCreate() 1759f7efa3c7SJed Brown @*/ 1760f7efa3c7SJed Brown PetscErrorCode ISLocalToGlobalMappingConcatenate(MPI_Comm comm,PetscInt n,const ISLocalToGlobalMapping ltogs[],ISLocalToGlobalMapping *ltogcat) 1761f7efa3c7SJed Brown { 1762f7efa3c7SJed Brown PetscInt i,cnt,m,*idx; 1763f7efa3c7SJed Brown PetscErrorCode ierr; 1764f7efa3c7SJed Brown 1765f7efa3c7SJed Brown PetscFunctionBegin; 1766f7efa3c7SJed Brown if (n < 0) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have a non-negative number of mappings, given %D",n); 1767f7efa3c7SJed Brown if (n > 0) PetscValidPointer(ltogs,3); 1768f7efa3c7SJed Brown for (i=0; i<n; i++) PetscValidHeaderSpecific(ltogs[i],IS_LTOGM_CLASSID,3); 1769f7efa3c7SJed Brown PetscValidPointer(ltogcat,4); 1770f7efa3c7SJed Brown for (cnt=0,i=0; i<n; i++) { 1771f7efa3c7SJed Brown ierr = ISLocalToGlobalMappingGetSize(ltogs[i],&m);CHKERRQ(ierr); 1772f7efa3c7SJed Brown cnt += m; 1773f7efa3c7SJed Brown } 1774785e854fSJed Brown ierr = PetscMalloc1(cnt,&idx);CHKERRQ(ierr); 1775f7efa3c7SJed Brown for (cnt=0,i=0; i<n; i++) { 1776f7efa3c7SJed Brown const PetscInt *subidx; 1777f7efa3c7SJed Brown ierr = ISLocalToGlobalMappingGetSize(ltogs[i],&m);CHKERRQ(ierr); 1778f7efa3c7SJed Brown ierr = ISLocalToGlobalMappingGetIndices(ltogs[i],&subidx);CHKERRQ(ierr); 1779f7efa3c7SJed Brown ierr = PetscMemcpy(&idx[cnt],subidx,m*sizeof(PetscInt));CHKERRQ(ierr); 1780f7efa3c7SJed Brown ierr = ISLocalToGlobalMappingRestoreIndices(ltogs[i],&subidx);CHKERRQ(ierr); 1781f7efa3c7SJed Brown cnt += m; 1782f7efa3c7SJed Brown } 1783f0413b6fSBarry Smith ierr = ISLocalToGlobalMappingCreate(comm,1,cnt,idx,PETSC_OWN_POINTER,ltogcat);CHKERRQ(ierr); 1784f7efa3c7SJed Brown PetscFunctionReturn(0); 1785f7efa3c7SJed Brown } 178604a59952SBarry Smith 1787413f72f0SBarry Smith /*MC 1788413f72f0SBarry Smith ISLOCALTOGLOBALMAPPINGBASIC - basic implementation of the ISLocalToGlobalMapping object. When ISGlobalToLocalMappingApply() is 1789413f72f0SBarry Smith used this is good for only small and moderate size problems. 1790413f72f0SBarry Smith 1791413f72f0SBarry Smith Options Database Keys: 1792413f72f0SBarry Smith + -islocaltoglobalmapping_type basic - select this method 1793413f72f0SBarry Smith 1794413f72f0SBarry Smith Level: beginner 1795413f72f0SBarry Smith 1796413f72f0SBarry Smith .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingSetType(), ISLOCALTOGLOBALMAPPINGHASH 1797413f72f0SBarry Smith M*/ 1798413f72f0SBarry Smith PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate_Basic(ISLocalToGlobalMapping ltog) 1799413f72f0SBarry Smith { 1800413f72f0SBarry Smith PetscFunctionBegin; 1801413f72f0SBarry Smith ltog->ops->globaltolocalmappingapply = ISGlobalToLocalMappingApply_Basic; 1802413f72f0SBarry Smith ltog->ops->globaltolocalmappingsetup = ISGlobalToLocalMappingSetUp_Basic; 1803413f72f0SBarry Smith ltog->ops->globaltolocalmappingapplyblock = ISGlobalToLocalMappingApplyBlock_Basic; 1804413f72f0SBarry Smith ltog->ops->destroy = ISLocalToGlobalMappingDestroy_Basic; 1805413f72f0SBarry Smith PetscFunctionReturn(0); 1806413f72f0SBarry Smith } 1807413f72f0SBarry Smith 1808413f72f0SBarry Smith /*MC 1809413f72f0SBarry Smith ISLOCALTOGLOBALMAPPINGHASH - hash implementation of the ISLocalToGlobalMapping object. When ISGlobalToLocalMappingApply() is 1810ed56e8eaSBarry Smith used this is good for large memory problems. 1811413f72f0SBarry Smith 1812413f72f0SBarry Smith Options Database Keys: 1813413f72f0SBarry Smith + -islocaltoglobalmapping_type hash - select this method 1814413f72f0SBarry Smith 1815ed56e8eaSBarry Smith 181695452b02SPatrick Sanan Notes: 181795452b02SPatrick Sanan This is selected automatically for large problems if the user does not set the type. 1818ed56e8eaSBarry Smith 1819413f72f0SBarry Smith Level: beginner 1820413f72f0SBarry Smith 1821413f72f0SBarry Smith .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingSetType(), ISLOCALTOGLOBALMAPPINGHASH 1822413f72f0SBarry Smith M*/ 1823413f72f0SBarry Smith PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate_Hash(ISLocalToGlobalMapping ltog) 1824413f72f0SBarry Smith { 1825413f72f0SBarry Smith PetscFunctionBegin; 1826413f72f0SBarry Smith ltog->ops->globaltolocalmappingapply = ISGlobalToLocalMappingApply_Hash; 1827413f72f0SBarry Smith ltog->ops->globaltolocalmappingsetup = ISGlobalToLocalMappingSetUp_Hash; 1828413f72f0SBarry Smith ltog->ops->globaltolocalmappingapplyblock = ISGlobalToLocalMappingApplyBlock_Hash; 1829413f72f0SBarry Smith ltog->ops->destroy = ISLocalToGlobalMappingDestroy_Hash; 1830413f72f0SBarry Smith PetscFunctionReturn(0); 1831413f72f0SBarry Smith } 1832413f72f0SBarry Smith 1833413f72f0SBarry Smith 1834413f72f0SBarry Smith /*@C 1835413f72f0SBarry Smith ISLocalToGlobalMappingRegister - Adds a method for applying a global to local mapping with an ISLocalToGlobalMapping 1836413f72f0SBarry Smith 1837413f72f0SBarry Smith Not Collective 1838413f72f0SBarry Smith 1839413f72f0SBarry Smith Input Parameters: 1840413f72f0SBarry Smith + sname - name of a new method 1841413f72f0SBarry Smith - routine_create - routine to create method context 1842413f72f0SBarry Smith 1843413f72f0SBarry Smith Notes: 1844ed56e8eaSBarry Smith ISLocalToGlobalMappingRegister() may be called multiple times to add several user-defined mappings. 1845413f72f0SBarry Smith 1846413f72f0SBarry Smith Sample usage: 1847413f72f0SBarry Smith .vb 1848413f72f0SBarry Smith ISLocalToGlobalMappingRegister("my_mapper",MyCreate); 1849413f72f0SBarry Smith .ve 1850413f72f0SBarry Smith 1851ed56e8eaSBarry Smith Then, your mapping can be chosen with the procedural interface via 1852413f72f0SBarry Smith $ ISLocalToGlobalMappingSetType(ltog,"my_mapper") 1853413f72f0SBarry Smith or at runtime via the option 1854ed56e8eaSBarry Smith $ -islocaltoglobalmapping_type my_mapper 1855413f72f0SBarry Smith 1856413f72f0SBarry Smith Level: advanced 1857413f72f0SBarry Smith 1858413f72f0SBarry Smith .keywords: ISLocalToGlobalMappingType, register 1859413f72f0SBarry Smith 1860413f72f0SBarry Smith .seealso: ISLocalToGlobalMappingRegisterAll(), ISLocalToGlobalMappingRegisterDestroy(), ISLOCALTOGLOBALMAPPINGBASIC, ISLOCALTOGLOBALMAPPINGHASH 1861413f72f0SBarry Smith 1862413f72f0SBarry Smith @*/ 1863413f72f0SBarry Smith PetscErrorCode ISLocalToGlobalMappingRegister(const char sname[],PetscErrorCode (*function)(ISLocalToGlobalMapping)) 1864413f72f0SBarry Smith { 1865413f72f0SBarry Smith PetscErrorCode ierr; 1866413f72f0SBarry Smith 1867413f72f0SBarry Smith PetscFunctionBegin; 1868413f72f0SBarry Smith ierr = PetscFunctionListAdd(&ISLocalToGlobalMappingList,sname,function);CHKERRQ(ierr); 1869413f72f0SBarry Smith PetscFunctionReturn(0); 1870413f72f0SBarry Smith } 1871413f72f0SBarry Smith 1872413f72f0SBarry Smith /*@C 1873ed56e8eaSBarry Smith ISLocalToGlobalMappingSetType - Builds ISLocalToGlobalMapping for a particular global to local mapping approach. 1874413f72f0SBarry Smith 1875413f72f0SBarry Smith Logically Collective on ISLocalToGlobalMapping 1876413f72f0SBarry Smith 1877413f72f0SBarry Smith Input Parameters: 1878413f72f0SBarry Smith + ltog - the ISLocalToGlobalMapping object 1879413f72f0SBarry Smith - type - a known method 1880413f72f0SBarry Smith 1881413f72f0SBarry Smith Options Database Key: 1882ed56e8eaSBarry Smith . -islocaltoglobalmapping_type <method> - Sets the method; use -help for a list 1883413f72f0SBarry Smith of available methods (for instance, basic or hash) 1884413f72f0SBarry Smith 1885413f72f0SBarry Smith Notes: 1886413f72f0SBarry Smith See "petsc/include/petscis.h" for available methods 1887413f72f0SBarry Smith 1888413f72f0SBarry Smith Normally, it is best to use the ISLocalToGlobalMappingSetFromOptions() command and 1889413f72f0SBarry Smith then set the ISLocalToGlobalMapping type from the options database rather than by using 1890413f72f0SBarry Smith this routine. 1891413f72f0SBarry Smith 1892413f72f0SBarry Smith Level: intermediate 1893413f72f0SBarry Smith 1894413f72f0SBarry Smith Developer Note: ISLocalToGlobalMappingRegister() is used to add new types to ISLocalToGlobalMappingList from which they 1895413f72f0SBarry Smith are accessed by ISLocalToGlobalMappingSetType(). 1896413f72f0SBarry Smith 1897413f72f0SBarry Smith .keywords: ISLocalToGlobalMapping, set, method 1898413f72f0SBarry Smith 1899413f72f0SBarry Smith .seealso: PCSetType(), ISLocalToGlobalMappingType, ISLocalToGlobalMappingRegister(), ISLocalToGlobalMappingCreate() 1900413f72f0SBarry Smith 1901413f72f0SBarry Smith @*/ 1902413f72f0SBarry Smith PetscErrorCode ISLocalToGlobalMappingSetType(ISLocalToGlobalMapping ltog, ISLocalToGlobalMappingType type) 1903413f72f0SBarry Smith { 1904413f72f0SBarry Smith PetscErrorCode ierr,(*r)(ISLocalToGlobalMapping); 1905413f72f0SBarry Smith PetscBool match; 1906413f72f0SBarry Smith 1907413f72f0SBarry Smith PetscFunctionBegin; 1908413f72f0SBarry Smith PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1); 1909413f72f0SBarry Smith PetscValidCharPointer(type,2); 1910413f72f0SBarry Smith 1911413f72f0SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)ltog,type,&match);CHKERRQ(ierr); 1912413f72f0SBarry Smith if (match) PetscFunctionReturn(0); 1913413f72f0SBarry Smith 1914413f72f0SBarry Smith ierr = PetscFunctionListFind(ISLocalToGlobalMappingList,type,&r);CHKERRQ(ierr); 1915413f72f0SBarry Smith if (!r) SETERRQ1(PetscObjectComm((PetscObject)ltog),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested ISLocalToGlobalMapping type %s",type); 1916413f72f0SBarry Smith /* Destroy the previous private LTOG context */ 1917413f72f0SBarry Smith if (ltog->ops->destroy) { 1918413f72f0SBarry Smith ierr = (*ltog->ops->destroy)(ltog);CHKERRQ(ierr); 1919413f72f0SBarry Smith ltog->ops->destroy = NULL; 1920413f72f0SBarry Smith } 1921413f72f0SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)ltog,type);CHKERRQ(ierr); 1922413f72f0SBarry Smith ierr = (*r)(ltog);CHKERRQ(ierr); 1923413f72f0SBarry Smith PetscFunctionReturn(0); 1924413f72f0SBarry Smith } 1925413f72f0SBarry Smith 1926413f72f0SBarry Smith PetscBool ISLocalToGlobalMappingRegisterAllCalled = PETSC_FALSE; 1927413f72f0SBarry Smith 1928413f72f0SBarry Smith /*@C 1929413f72f0SBarry Smith ISLocalToGlobalMappingRegisterAll - Registers all of the local to global mapping components in the IS package. 1930413f72f0SBarry Smith 1931413f72f0SBarry Smith Not Collective 1932413f72f0SBarry Smith 1933413f72f0SBarry Smith Level: advanced 1934413f72f0SBarry Smith 1935413f72f0SBarry Smith .keywords: IS, register, all 1936413f72f0SBarry Smith .seealso: ISRegister(), ISLocalToGlobalRegister() 1937413f72f0SBarry Smith @*/ 1938413f72f0SBarry Smith PetscErrorCode ISLocalToGlobalMappingRegisterAll(void) 1939413f72f0SBarry Smith { 1940413f72f0SBarry Smith PetscErrorCode ierr; 1941413f72f0SBarry Smith 1942413f72f0SBarry Smith PetscFunctionBegin; 1943413f72f0SBarry Smith if (ISLocalToGlobalMappingRegisterAllCalled) PetscFunctionReturn(0); 1944413f72f0SBarry Smith ISLocalToGlobalMappingRegisterAllCalled = PETSC_TRUE; 1945413f72f0SBarry Smith 1946413f72f0SBarry Smith ierr = ISLocalToGlobalMappingRegister(ISLOCALTOGLOBALMAPPINGBASIC, ISLocalToGlobalMappingCreate_Basic);CHKERRQ(ierr); 1947413f72f0SBarry Smith ierr = ISLocalToGlobalMappingRegister(ISLOCALTOGLOBALMAPPINGHASH, ISLocalToGlobalMappingCreate_Hash);CHKERRQ(ierr); 1948413f72f0SBarry Smith PetscFunctionReturn(0); 1949413f72f0SBarry Smith } 195004a59952SBarry Smith 1951