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 199305a4c7SMatthew G. Knepley PetscErrorCode ISGetPointRange(IS pointIS, PetscInt *pStart, PetscInt *pEnd, const PetscInt **points) 209305a4c7SMatthew G. Knepley { 219305a4c7SMatthew G. Knepley PetscInt numCells, step = 1; 229305a4c7SMatthew G. Knepley PetscBool isStride; 239305a4c7SMatthew G. Knepley PetscErrorCode ierr; 249305a4c7SMatthew G. Knepley 259305a4c7SMatthew G. Knepley PetscFunctionBeginHot; 269305a4c7SMatthew G. Knepley *pStart = 0; 279305a4c7SMatthew G. Knepley *points = NULL; 289305a4c7SMatthew G. Knepley ierr = ISGetLocalSize(pointIS, &numCells);CHKERRQ(ierr); 299305a4c7SMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject) pointIS, ISSTRIDE, &isStride);CHKERRQ(ierr); 309305a4c7SMatthew G. Knepley if (isStride) {ierr = ISStrideGetInfo(pointIS, pStart, &step);CHKERRQ(ierr);} 319305a4c7SMatthew G. Knepley *pEnd = *pStart + numCells; 329305a4c7SMatthew G. Knepley if (!isStride || step != 1) {ierr = ISGetIndices(pointIS, points);CHKERRQ(ierr);} 339305a4c7SMatthew G. Knepley PetscFunctionReturn(0); 349305a4c7SMatthew G. Knepley } 359305a4c7SMatthew G. Knepley 369305a4c7SMatthew G. Knepley PetscErrorCode ISRestorePointRange(IS pointIS, PetscInt *pStart, PetscInt *pEnd, const PetscInt **points) 379305a4c7SMatthew G. Knepley { 389305a4c7SMatthew G. Knepley PetscInt step = 1; 399305a4c7SMatthew G. Knepley PetscBool isStride; 409305a4c7SMatthew G. Knepley PetscErrorCode ierr; 419305a4c7SMatthew G. Knepley 429305a4c7SMatthew G. Knepley PetscFunctionBeginHot; 439305a4c7SMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject) pointIS, ISSTRIDE, &isStride);CHKERRQ(ierr); 449305a4c7SMatthew G. Knepley if (isStride) {ierr = ISStrideGetInfo(pointIS, pStart, &step);CHKERRQ(ierr);} 459305a4c7SMatthew G. Knepley if (!isStride || step != 1) {ierr = ISGetIndices(pointIS, points);CHKERRQ(ierr);} 469305a4c7SMatthew G. Knepley PetscFunctionReturn(0); 479305a4c7SMatthew G. Knepley } 489305a4c7SMatthew G. Knepley 499305a4c7SMatthew G. Knepley PetscErrorCode ISGetPointSubrange(IS subpointIS, PetscInt pStart, PetscInt pEnd, const PetscInt *points) 509305a4c7SMatthew G. Knepley { 519305a4c7SMatthew G. Knepley PetscErrorCode ierr; 529305a4c7SMatthew G. Knepley 539305a4c7SMatthew G. Knepley PetscFunctionBeginHot; 549305a4c7SMatthew G. Knepley if (points) { 559305a4c7SMatthew G. Knepley ierr = ISSetType(subpointIS, ISGENERAL);CHKERRQ(ierr); 569305a4c7SMatthew G. Knepley ierr = ISGeneralSetIndices(subpointIS, pEnd-pStart, &points[pStart], PETSC_USE_POINTER);CHKERRQ(ierr); 579305a4c7SMatthew G. Knepley } else { 589305a4c7SMatthew G. Knepley ierr = ISSetType(subpointIS, ISSTRIDE);CHKERRQ(ierr); 599305a4c7SMatthew G. Knepley ierr = ISStrideSetStride(subpointIS, pEnd-pStart, pStart, 1);CHKERRQ(ierr); 609305a4c7SMatthew G. Knepley } 619305a4c7SMatthew G. Knepley PetscFunctionReturn(0); 629305a4c7SMatthew G. Knepley } 639305a4c7SMatthew 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: 396*a2b725a8SWilliam Gropp + mapping - mapping data structure 397*a2b725a8SWilliam Gropp - 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 6474cb36875SStefano Zampini Collective on is 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 6564cb36875SStefano Zampini Notes: 6574cb36875SStefano Zampini The output IS will have the same communicator of the input IS. 6584cb36875SStefano Zampini 659a997ad1aSLois Curfman McInnes Level: advanced 660a997ad1aSLois Curfman McInnes 661273d9f13SBarry Smith Concepts: mapping^local to global 6623acfe500SLois Curfman McInnes 66390f02eecSBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingCreate(), 664d4bb536fSBarry Smith ISLocalToGlobalMappingDestroy(), ISGlobalToLocalMappingApply() 66590f02eecSBarry Smith @*/ 6667087cfbeSBarry Smith PetscErrorCode ISLocalToGlobalMappingApplyIS(ISLocalToGlobalMapping mapping,IS is,IS *newis) 66790f02eecSBarry Smith { 6686849ba73SBarry Smith PetscErrorCode ierr; 669e24637baSBarry Smith PetscInt n,*idxout; 6705d0c19d7SBarry Smith const PetscInt *idxin; 6713a40ed3dSBarry Smith 6723a40ed3dSBarry Smith PetscFunctionBegin; 6730700a824SBarry Smith PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1); 6740700a824SBarry Smith PetscValidHeaderSpecific(is,IS_CLASSID,2); 6754482741eSBarry Smith PetscValidPointer(newis,3); 67690f02eecSBarry Smith 6773b9aefa3SBarry Smith ierr = ISGetLocalSize(is,&n);CHKERRQ(ierr); 67890f02eecSBarry Smith ierr = ISGetIndices(is,&idxin);CHKERRQ(ierr); 679785e854fSJed Brown ierr = PetscMalloc1(n,&idxout);CHKERRQ(ierr); 680e24637baSBarry Smith ierr = ISLocalToGlobalMappingApply(mapping,n,idxin,idxout);CHKERRQ(ierr); 6813b9aefa3SBarry Smith ierr = ISRestoreIndices(is,&idxin);CHKERRQ(ierr); 682543f3098SMatthew G. Knepley ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is),n,idxout,PETSC_OWN_POINTER,newis);CHKERRQ(ierr); 6833a40ed3dSBarry Smith PetscFunctionReturn(0); 68490f02eecSBarry Smith } 68590f02eecSBarry Smith 686b89cb25eSSatish Balay /*@ 6873acfe500SLois Curfman McInnes ISLocalToGlobalMappingApply - Takes a list of integers in a local numbering 6883acfe500SLois Curfman McInnes and converts them to the global numbering. 68990f02eecSBarry Smith 690b9cd556bSLois Curfman McInnes Not collective 691b9cd556bSLois Curfman McInnes 692bb25748dSBarry Smith Input Parameters: 693b9cd556bSLois Curfman McInnes + mapping - the local to global mapping context 694bb25748dSBarry Smith . N - number of integers 695b9cd556bSLois Curfman McInnes - in - input indices in local numbering 696bb25748dSBarry Smith 697bb25748dSBarry Smith Output Parameter: 698bb25748dSBarry Smith . out - indices in global numbering 699bb25748dSBarry Smith 700b9cd556bSLois Curfman McInnes Notes: 701b9cd556bSLois Curfman McInnes The in and out array parameters may be identical. 702d4bb536fSBarry Smith 703a997ad1aSLois Curfman McInnes Level: advanced 704a997ad1aSLois Curfman McInnes 70545b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingApplyBlock(), ISLocalToGlobalMappingCreate(),ISLocalToGlobalMappingDestroy(), 7060752156aSBarry Smith ISLocalToGlobalMappingApplyIS(),AOCreateBasic(),AOApplicationToPetsc(), 707d4bb536fSBarry Smith AOPetscToApplication(), ISGlobalToLocalMappingApply() 708bb25748dSBarry Smith 709273d9f13SBarry Smith Concepts: mapping^local to global 710afcb2eb5SJed Brown @*/ 711afcb2eb5SJed Brown PetscErrorCode ISLocalToGlobalMappingApply(ISLocalToGlobalMapping mapping,PetscInt N,const PetscInt in[],PetscInt out[]) 712afcb2eb5SJed Brown { 713cbc1caf0SMatthew G. Knepley PetscInt i,bs,Nmax; 71445b6f7e9SBarry Smith 71545b6f7e9SBarry Smith PetscFunctionBegin; 716cbc1caf0SMatthew G. Knepley PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1); 717cbc1caf0SMatthew G. Knepley bs = mapping->bs; 718cbc1caf0SMatthew G. Knepley Nmax = bs*mapping->n; 71945b6f7e9SBarry Smith if (bs == 1) { 720cbc1caf0SMatthew G. Knepley const PetscInt *idx = mapping->indices; 72145b6f7e9SBarry Smith for (i=0; i<N; i++) { 72245b6f7e9SBarry Smith if (in[i] < 0) { 72345b6f7e9SBarry Smith out[i] = in[i]; 72445b6f7e9SBarry Smith continue; 72545b6f7e9SBarry Smith } 726e24637baSBarry 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); 72745b6f7e9SBarry Smith out[i] = idx[in[i]]; 72845b6f7e9SBarry Smith } 72945b6f7e9SBarry Smith } else { 730cbc1caf0SMatthew G. Knepley const PetscInt *idx = mapping->indices; 73145b6f7e9SBarry Smith for (i=0; i<N; i++) { 73245b6f7e9SBarry Smith if (in[i] < 0) { 73345b6f7e9SBarry Smith out[i] = in[i]; 73445b6f7e9SBarry Smith continue; 73545b6f7e9SBarry Smith } 736e24637baSBarry 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); 73745b6f7e9SBarry Smith out[i] = idx[in[i]/bs]*bs + (in[i] % bs); 73845b6f7e9SBarry Smith } 73945b6f7e9SBarry Smith } 74045b6f7e9SBarry Smith PetscFunctionReturn(0); 74145b6f7e9SBarry Smith } 74245b6f7e9SBarry Smith 74345b6f7e9SBarry Smith /*@ 7446006e8d2SBarry Smith ISLocalToGlobalMappingApplyBlock - Takes a list of integers in a local block numbering and converts them to the global block numbering 74545b6f7e9SBarry Smith 74645b6f7e9SBarry Smith Not collective 74745b6f7e9SBarry Smith 74845b6f7e9SBarry Smith Input Parameters: 74945b6f7e9SBarry Smith + mapping - the local to global mapping context 75045b6f7e9SBarry Smith . N - number of integers 7516006e8d2SBarry Smith - in - input indices in local block numbering 75245b6f7e9SBarry Smith 75345b6f7e9SBarry Smith Output Parameter: 7546006e8d2SBarry Smith . out - indices in global block numbering 75545b6f7e9SBarry Smith 75645b6f7e9SBarry Smith Notes: 75745b6f7e9SBarry Smith The in and out array parameters may be identical. 75845b6f7e9SBarry Smith 7596006e8d2SBarry Smith Example: 7606006e8d2SBarry 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 7616006e8d2SBarry Smith (the first block) would produce 0 and the mapping applied to 1 (the second block) would produce 3. 7626006e8d2SBarry Smith 76345b6f7e9SBarry Smith Level: advanced 76445b6f7e9SBarry Smith 76545b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingCreate(),ISLocalToGlobalMappingDestroy(), 76645b6f7e9SBarry Smith ISLocalToGlobalMappingApplyIS(),AOCreateBasic(),AOApplicationToPetsc(), 76745b6f7e9SBarry Smith AOPetscToApplication(), ISGlobalToLocalMappingApply() 76845b6f7e9SBarry Smith 76945b6f7e9SBarry Smith Concepts: mapping^local to global 77045b6f7e9SBarry Smith @*/ 77145b6f7e9SBarry Smith PetscErrorCode ISLocalToGlobalMappingApplyBlock(ISLocalToGlobalMapping mapping,PetscInt N,const PetscInt in[],PetscInt out[]) 77245b6f7e9SBarry Smith { 773cbc1caf0SMatthew G. Knepley 774cbc1caf0SMatthew G. Knepley PetscFunctionBegin; 775cbc1caf0SMatthew G. Knepley PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1); 776cbc1caf0SMatthew G. Knepley { 777afcb2eb5SJed Brown PetscInt i,Nmax = mapping->n; 778afcb2eb5SJed Brown const PetscInt *idx = mapping->indices; 779d4bb536fSBarry Smith 780afcb2eb5SJed Brown for (i=0; i<N; i++) { 781afcb2eb5SJed Brown if (in[i] < 0) { 782afcb2eb5SJed Brown out[i] = in[i]; 783afcb2eb5SJed Brown continue; 784afcb2eb5SJed Brown } 785e24637baSBarry 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); 786afcb2eb5SJed Brown out[i] = idx[in[i]]; 787afcb2eb5SJed Brown } 788cbc1caf0SMatthew G. Knepley } 789afcb2eb5SJed Brown PetscFunctionReturn(0); 790afcb2eb5SJed Brown } 791d4bb536fSBarry Smith 7927e99dc12SLawrence Mitchell /*@ 793a997ad1aSLois Curfman McInnes ISGlobalToLocalMappingApply - Provides the local numbering for a list of integers 794a997ad1aSLois Curfman McInnes specified with a global numbering. 795d4bb536fSBarry Smith 796b9cd556bSLois Curfman McInnes Not collective 797b9cd556bSLois Curfman McInnes 798d4bb536fSBarry Smith Input Parameters: 799b9cd556bSLois Curfman McInnes + mapping - mapping between local and global numbering 800d4bb536fSBarry Smith . type - IS_GTOLM_MASK - replaces global indices with no local value with -1 801d4bb536fSBarry Smith IS_GTOLM_DROP - drops the indices with no local value from the output list 802d4bb536fSBarry Smith . n - number of global indices to map 803b9cd556bSLois Curfman McInnes - idx - global indices to map 804d4bb536fSBarry Smith 805d4bb536fSBarry Smith Output Parameters: 806b9cd556bSLois Curfman McInnes + nout - number of indices in output array (if type == IS_GTOLM_MASK then nout = n) 807b9cd556bSLois Curfman McInnes - idxout - local index of each global index, one must pass in an array long enough 808e182c471SBarry Smith to hold all the indices. You can call ISGlobalToLocalMappingApply() with 8090298fd71SBarry Smith idxout == NULL to determine the required length (returned in nout) 810e182c471SBarry Smith and then allocate the required space and call ISGlobalToLocalMappingApply() 811e182c471SBarry Smith a second time to set the values. 812d4bb536fSBarry Smith 813b9cd556bSLois Curfman McInnes Notes: 8140298fd71SBarry Smith Either nout or idxout may be NULL. idx and idxout may be identical. 815d4bb536fSBarry Smith 8169a7b7924SJed Brown For "small" problems when using ISGlobalToLocalMappingApply() and ISGlobalToLocalMappingApplyBlock(), the ISLocalToGlobalMappingType of ISLOCALTOGLOBALMAPPINGBASIC will be used; 817413f72f0SBarry 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. 818413f72f0SBarry Smith Use ISLocalToGlobalMappingSetType() or call ISLocalToGlobalMappingSetFromOptions() with the option -islocaltoglobalmapping_type <basic,hash> to control which is used. 8190f5bd95cSBarry Smith 820a997ad1aSLois Curfman McInnes Level: advanced 821a997ad1aSLois Curfman McInnes 82232fd6b96SBarry Smith Developer Note: The manual page states that idx and idxout may be identical but the calling 82332fd6b96SBarry Smith sequence declares idx as const so it cannot be the same as idxout. 82432fd6b96SBarry Smith 825273d9f13SBarry Smith Concepts: mapping^global to local 826d4bb536fSBarry Smith 8279d90f715SBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISGlobalToLocalMappingApplyBlock(), ISLocalToGlobalMappingCreate(), 828413f72f0SBarry Smith ISLocalToGlobalMappingDestroy() 829d4bb536fSBarry Smith @*/ 830413f72f0SBarry Smith PetscErrorCode ISGlobalToLocalMappingApply(ISLocalToGlobalMapping mapping,ISGlobalToLocalMappingMode type,PetscInt n,const PetscInt idx[],PetscInt *nout,PetscInt idxout[]) 831d4bb536fSBarry Smith { 8329d90f715SBarry Smith PetscErrorCode ierr; 8339d90f715SBarry Smith 8349d90f715SBarry Smith PetscFunctionBegin; 8359d90f715SBarry Smith PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1); 836413f72f0SBarry Smith if (!mapping->data) { 837413f72f0SBarry Smith ierr = ISGlobalToLocalMappingSetUp(mapping);CHKERRQ(ierr); 8389d90f715SBarry Smith } 839413f72f0SBarry Smith ierr = (*mapping->ops->globaltolocalmappingapply)(mapping,type,n,idx,nout,idxout);CHKERRQ(ierr); 8409d90f715SBarry Smith PetscFunctionReturn(0); 8419d90f715SBarry Smith } 8429d90f715SBarry Smith 843d4fe737eSStefano Zampini /*@ 844d4fe737eSStefano Zampini ISGlobalToLocalMappingApplyIS - Creates from an IS in the global numbering 845d4fe737eSStefano Zampini a new index set using the local numbering defined in an ISLocalToGlobalMapping 846d4fe737eSStefano Zampini context. 847d4fe737eSStefano Zampini 848d4fe737eSStefano Zampini Not collective 849d4fe737eSStefano Zampini 850d4fe737eSStefano Zampini Input Parameters: 851d4fe737eSStefano Zampini + mapping - mapping between local and global numbering 8522785b321SStefano Zampini . type - IS_GTOLM_MASK - replaces global indices with no local value with -1 8532785b321SStefano Zampini IS_GTOLM_DROP - drops the indices with no local value from the output list 854d4fe737eSStefano Zampini - is - index set in global numbering 855d4fe737eSStefano Zampini 856d4fe737eSStefano Zampini Output Parameters: 857d4fe737eSStefano Zampini . newis - index set in local numbering 858d4fe737eSStefano Zampini 8594cb36875SStefano Zampini Notes: 8604cb36875SStefano Zampini The output IS will be sequential, as it encodes a purely local operation 8614cb36875SStefano Zampini 862d4fe737eSStefano Zampini Level: advanced 863d4fe737eSStefano Zampini 864d4fe737eSStefano Zampini Concepts: mapping^local to global 865d4fe737eSStefano Zampini 866d4fe737eSStefano Zampini .seealso: ISGlobalToLocalMappingApply(), ISLocalToGlobalMappingCreate(), 867d4fe737eSStefano Zampini ISLocalToGlobalMappingDestroy() 868d4fe737eSStefano Zampini @*/ 869413f72f0SBarry Smith PetscErrorCode ISGlobalToLocalMappingApplyIS(ISLocalToGlobalMapping mapping,ISGlobalToLocalMappingMode type,IS is,IS *newis) 870d4fe737eSStefano Zampini { 871d4fe737eSStefano Zampini PetscErrorCode ierr; 872d4fe737eSStefano Zampini PetscInt n,nout,*idxout; 873d4fe737eSStefano Zampini const PetscInt *idxin; 874d4fe737eSStefano Zampini 875d4fe737eSStefano Zampini PetscFunctionBegin; 876d4fe737eSStefano Zampini PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1); 877d4fe737eSStefano Zampini PetscValidHeaderSpecific(is,IS_CLASSID,3); 878d4fe737eSStefano Zampini PetscValidPointer(newis,4); 879d4fe737eSStefano Zampini 880d4fe737eSStefano Zampini ierr = ISGetLocalSize(is,&n);CHKERRQ(ierr); 881d4fe737eSStefano Zampini ierr = ISGetIndices(is,&idxin);CHKERRQ(ierr); 882d4fe737eSStefano Zampini if (type == IS_GTOLM_MASK) { 883d4fe737eSStefano Zampini ierr = PetscMalloc1(n,&idxout);CHKERRQ(ierr); 884d4fe737eSStefano Zampini } else { 885d4fe737eSStefano Zampini ierr = ISGlobalToLocalMappingApply(mapping,type,n,idxin,&nout,NULL);CHKERRQ(ierr); 886d4fe737eSStefano Zampini ierr = PetscMalloc1(nout,&idxout);CHKERRQ(ierr); 887d4fe737eSStefano Zampini } 888d4fe737eSStefano Zampini ierr = ISGlobalToLocalMappingApply(mapping,type,n,idxin,&nout,idxout);CHKERRQ(ierr); 889d4fe737eSStefano Zampini ierr = ISRestoreIndices(is,&idxin);CHKERRQ(ierr); 8904cb36875SStefano Zampini ierr = ISCreateGeneral(PETSC_COMM_SELF,nout,idxout,PETSC_OWN_POINTER,newis);CHKERRQ(ierr); 891d4fe737eSStefano Zampini PetscFunctionReturn(0); 892d4fe737eSStefano Zampini } 893d4fe737eSStefano Zampini 8949d90f715SBarry Smith /*@ 8959d90f715SBarry Smith ISGlobalToLocalMappingApplyBlock - Provides the local block numbering for a list of integers 8969d90f715SBarry Smith specified with a block global numbering. 8979d90f715SBarry Smith 8989d90f715SBarry Smith Not collective 8999d90f715SBarry Smith 9009d90f715SBarry Smith Input Parameters: 9019d90f715SBarry Smith + mapping - mapping between local and global numbering 9029d90f715SBarry Smith . type - IS_GTOLM_MASK - replaces global indices with no local value with -1 9039d90f715SBarry Smith IS_GTOLM_DROP - drops the indices with no local value from the output list 9049d90f715SBarry Smith . n - number of global indices to map 9059d90f715SBarry Smith - idx - global indices to map 9069d90f715SBarry Smith 9079d90f715SBarry Smith Output Parameters: 9089d90f715SBarry Smith + nout - number of indices in output array (if type == IS_GTOLM_MASK then nout = n) 9099d90f715SBarry Smith - idxout - local index of each global index, one must pass in an array long enough 9109d90f715SBarry Smith to hold all the indices. You can call ISGlobalToLocalMappingApplyBlock() with 9119d90f715SBarry Smith idxout == NULL to determine the required length (returned in nout) 9129d90f715SBarry Smith and then allocate the required space and call ISGlobalToLocalMappingApplyBlock() 9139d90f715SBarry Smith a second time to set the values. 9149d90f715SBarry Smith 9159d90f715SBarry Smith Notes: 9169d90f715SBarry Smith Either nout or idxout may be NULL. idx and idxout may be identical. 9179d90f715SBarry Smith 9189a7b7924SJed Brown For "small" problems when using ISGlobalToLocalMappingApply() and ISGlobalToLocalMappingApplyBlock(), the ISLocalToGlobalMappingType of ISLOCALTOGLOBALMAPPINGBASIC will be used; 919413f72f0SBarry 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. 920413f72f0SBarry Smith Use ISLocalToGlobalMappingSetType() or call ISLocalToGlobalMappingSetFromOptions() with the option -islocaltoglobalmapping_type <basic,hash> to control which is used. 9219d90f715SBarry Smith 9229d90f715SBarry Smith Level: advanced 9239d90f715SBarry Smith 9249d90f715SBarry Smith Developer Note: The manual page states that idx and idxout may be identical but the calling 9259d90f715SBarry Smith sequence declares idx as const so it cannot be the same as idxout. 9269d90f715SBarry Smith 9279d90f715SBarry Smith Concepts: mapping^global to local 9289d90f715SBarry Smith 9299d90f715SBarry Smith .seealso: ISLocalToGlobalMappingApply(), ISGlobalToLocalMappingApply(), ISLocalToGlobalMappingCreate(), 930413f72f0SBarry Smith ISLocalToGlobalMappingDestroy() 9319d90f715SBarry Smith @*/ 932413f72f0SBarry Smith PetscErrorCode ISGlobalToLocalMappingApplyBlock(ISLocalToGlobalMapping mapping,ISGlobalToLocalMappingMode type, 9339d90f715SBarry Smith PetscInt n,const PetscInt idx[],PetscInt *nout,PetscInt idxout[]) 9349d90f715SBarry Smith { 9356849ba73SBarry Smith PetscErrorCode ierr; 936d4bb536fSBarry Smith 9373a40ed3dSBarry Smith PetscFunctionBegin; 9380700a824SBarry Smith PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1); 939413f72f0SBarry Smith if (!mapping->data) { 940413f72f0SBarry Smith ierr = ISGlobalToLocalMappingSetUp(mapping);CHKERRQ(ierr); 941d4bb536fSBarry Smith } 942413f72f0SBarry Smith ierr = (*mapping->ops->globaltolocalmappingapplyblock)(mapping,type,n,idx,nout,idxout);CHKERRQ(ierr); 9433a40ed3dSBarry Smith PetscFunctionReturn(0); 944d4bb536fSBarry Smith } 94590f02eecSBarry Smith 946413f72f0SBarry Smith 94789d82c54SBarry Smith /*@C 9486a818285SBarry Smith ISLocalToGlobalMappingGetBlockInfo - Gets the neighbor information for each processor and 94989d82c54SBarry Smith each index shared by more than one processor 95089d82c54SBarry Smith 95189d82c54SBarry Smith Collective on ISLocalToGlobalMapping 95289d82c54SBarry Smith 95389d82c54SBarry Smith Input Parameters: 95489d82c54SBarry Smith . mapping - the mapping from local to global indexing 95589d82c54SBarry Smith 95689d82c54SBarry Smith Output Parameter: 95789d82c54SBarry Smith + nproc - number of processors that are connected to this one 95889d82c54SBarry Smith . proc - neighboring processors 95907b52d57SBarry Smith . numproc - number of indices for each subdomain (processor) 9603463a7baSJed Brown - indices - indices of nodes (in local numbering) shared with neighbors (sorted by global numbering) 96189d82c54SBarry Smith 96289d82c54SBarry Smith Level: advanced 96389d82c54SBarry Smith 964273d9f13SBarry Smith Concepts: mapping^local to global 96589d82c54SBarry Smith 9662cfcea29SBarry Smith Fortran Usage: 9672cfcea29SBarry Smith $ ISLocalToGlobalMpngGetInfoSize(ISLocalToGlobalMapping,PetscInt nproc,PetscInt numprocmax,ierr) followed by 9682cfcea29SBarry Smith $ ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt nproc, PetscInt procs[nproc],PetscInt numprocs[nproc], 9692cfcea29SBarry Smith PetscInt indices[nproc][numprocmax],ierr) 9702cfcea29SBarry Smith There is no ISLocalToGlobalMappingRestoreInfo() in Fortran. You must make sure that procs[], numprocs[] and 9712cfcea29SBarry Smith indices[][] are large enough arrays, either by allocating them dynamically or defining static ones large enough. 9722cfcea29SBarry Smith 9732cfcea29SBarry Smith 97407b52d57SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(), 97507b52d57SBarry Smith ISLocalToGlobalMappingRestoreInfo() 97689d82c54SBarry Smith @*/ 9776a818285SBarry Smith PetscErrorCode ISLocalToGlobalMappingGetBlockInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[]) 97889d82c54SBarry Smith { 9796849ba73SBarry Smith PetscErrorCode ierr; 980268a049cSStefano Zampini 981268a049cSStefano Zampini PetscFunctionBegin; 982268a049cSStefano Zampini PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1); 983268a049cSStefano Zampini if (mapping->info_cached) { 984268a049cSStefano Zampini *nproc = mapping->info_nproc; 985268a049cSStefano Zampini *procs = mapping->info_procs; 986268a049cSStefano Zampini *numprocs = mapping->info_numprocs; 987268a049cSStefano Zampini *indices = mapping->info_indices; 988268a049cSStefano Zampini } else { 989268a049cSStefano Zampini ierr = ISLocalToGlobalMappingGetBlockInfo_Private(mapping,nproc,procs,numprocs,indices);CHKERRQ(ierr); 990268a049cSStefano Zampini } 991268a049cSStefano Zampini PetscFunctionReturn(0); 992268a049cSStefano Zampini } 993268a049cSStefano Zampini 994268a049cSStefano Zampini static PetscErrorCode ISLocalToGlobalMappingGetBlockInfo_Private(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[]) 995268a049cSStefano Zampini { 996268a049cSStefano Zampini PetscErrorCode ierr; 99797f1f81fSBarry Smith PetscMPIInt size,rank,tag1,tag2,tag3,*len,*source,imdex; 99832dcc486SBarry Smith PetscInt i,n = mapping->n,Ng,ng,max = 0,*lindices = mapping->indices; 99932dcc486SBarry Smith PetscInt *nprocs,*owner,nsends,*sends,j,*starts,nmax,nrecvs,*recvs,proc; 100097f1f81fSBarry Smith PetscInt cnt,scale,*ownedsenders,*nownedsenders,rstart,nowned; 100132dcc486SBarry Smith PetscInt node,nownedm,nt,*sends2,nsends2,*starts2,*lens2,*dest,nrecvs2,*starts3,*recvs2,k,*bprocs,*tmp; 100232dcc486SBarry Smith PetscInt first_procs,first_numprocs,*first_indices; 100389d82c54SBarry Smith MPI_Request *recv_waits,*send_waits; 100430dcb7c9SBarry Smith MPI_Status recv_status,*send_status,*recv_statuses; 1005ce94432eSBarry Smith MPI_Comm comm; 1006ace3abfcSBarry Smith PetscBool debug = PETSC_FALSE; 100789d82c54SBarry Smith 100889d82c54SBarry Smith PetscFunctionBegin; 1009ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)mapping,&comm);CHKERRQ(ierr); 101024cf384cSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 101124cf384cSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 101224cf384cSBarry Smith if (size == 1) { 101324cf384cSBarry Smith *nproc = 0; 10140298fd71SBarry Smith *procs = NULL; 101595dccacaSBarry Smith ierr = PetscNew(numprocs);CHKERRQ(ierr); 10161e2105dcSBarry Smith (*numprocs)[0] = 0; 101795dccacaSBarry Smith ierr = PetscNew(indices);CHKERRQ(ierr); 10180298fd71SBarry Smith (*indices)[0] = NULL; 1019268a049cSStefano Zampini /* save info for reuse */ 1020268a049cSStefano Zampini mapping->info_nproc = *nproc; 1021268a049cSStefano Zampini mapping->info_procs = *procs; 1022268a049cSStefano Zampini mapping->info_numprocs = *numprocs; 1023268a049cSStefano Zampini mapping->info_indices = *indices; 1024268a049cSStefano Zampini mapping->info_cached = PETSC_TRUE; 102524cf384cSBarry Smith PetscFunctionReturn(0); 102624cf384cSBarry Smith } 102724cf384cSBarry Smith 1028c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)mapping)->options,NULL,"-islocaltoglobalmappinggetinfo_debug",&debug,NULL);CHKERRQ(ierr); 102907b52d57SBarry Smith 10303677ff5aSBarry Smith /* 10316a818285SBarry Smith Notes on ISLocalToGlobalMappingGetBlockInfo 10323677ff5aSBarry Smith 10333677ff5aSBarry Smith globally owned node - the nodes that have been assigned to this processor in global 10343677ff5aSBarry Smith numbering, just for this routine. 10353677ff5aSBarry Smith 10363677ff5aSBarry Smith nontrivial globally owned node - node assigned to this processor that is on a subdomain 10373677ff5aSBarry Smith boundary (i.e. is has more than one local owner) 10383677ff5aSBarry Smith 10393677ff5aSBarry Smith locally owned node - node that exists on this processors subdomain 10403677ff5aSBarry Smith 10413677ff5aSBarry Smith nontrivial locally owned node - node that is not in the interior (i.e. has more than one 10423677ff5aSBarry Smith local subdomain 10433677ff5aSBarry Smith */ 104424cf384cSBarry Smith ierr = PetscObjectGetNewTag((PetscObject)mapping,&tag1);CHKERRQ(ierr); 104524cf384cSBarry Smith ierr = PetscObjectGetNewTag((PetscObject)mapping,&tag2);CHKERRQ(ierr); 104624cf384cSBarry Smith ierr = PetscObjectGetNewTag((PetscObject)mapping,&tag3);CHKERRQ(ierr); 104789d82c54SBarry Smith 104889d82c54SBarry Smith for (i=0; i<n; i++) { 104989d82c54SBarry Smith if (lindices[i] > max) max = lindices[i]; 105089d82c54SBarry Smith } 1051b2566f29SBarry Smith ierr = MPIU_Allreduce(&max,&Ng,1,MPIU_INT,MPI_MAX,comm);CHKERRQ(ierr); 105278058e43SBarry Smith Ng++; 105389d82c54SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 105489d82c54SBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1055bc8ff85bSBarry Smith scale = Ng/size + 1; 1056a2e34c3dSBarry Smith ng = scale; if (rank == size-1) ng = Ng - scale*(size-1); ng = PetscMax(1,ng); 1057caba0dd0SBarry Smith rstart = scale*rank; 105889d82c54SBarry Smith 105989d82c54SBarry Smith /* determine ownership ranges of global indices */ 1060785e854fSJed Brown ierr = PetscMalloc1(2*size,&nprocs);CHKERRQ(ierr); 106132dcc486SBarry Smith ierr = PetscMemzero(nprocs,2*size*sizeof(PetscInt));CHKERRQ(ierr); 106289d82c54SBarry Smith 106389d82c54SBarry Smith /* determine owners of each local node */ 1064785e854fSJed Brown ierr = PetscMalloc1(n,&owner);CHKERRQ(ierr); 106589d82c54SBarry Smith for (i=0; i<n; i++) { 10663677ff5aSBarry Smith proc = lindices[i]/scale; /* processor that globally owns this index */ 106727c402fcSBarry Smith nprocs[2*proc+1] = 1; /* processor globally owns at least one of ours */ 10683677ff5aSBarry Smith owner[i] = proc; 106927c402fcSBarry Smith nprocs[2*proc]++; /* count of how many that processor globally owns of ours */ 107089d82c54SBarry Smith } 107127c402fcSBarry Smith nsends = 0; for (i=0; i<size; i++) nsends += nprocs[2*i+1]; 10727904a332SBarry Smith ierr = PetscInfo1(mapping,"Number of global owners for my local data %D\n",nsends);CHKERRQ(ierr); 107389d82c54SBarry Smith 107489d82c54SBarry Smith /* inform other processors of number of messages and max length*/ 107527c402fcSBarry Smith ierr = PetscMaxSum(comm,nprocs,&nmax,&nrecvs);CHKERRQ(ierr); 10767904a332SBarry Smith ierr = PetscInfo1(mapping,"Number of local owners for my global data %D\n",nrecvs);CHKERRQ(ierr); 107789d82c54SBarry Smith 107889d82c54SBarry Smith /* post receives for owned rows */ 1079785e854fSJed Brown ierr = PetscMalloc1((2*nrecvs+1)*(nmax+1),&recvs);CHKERRQ(ierr); 1080854ce69bSBarry Smith ierr = PetscMalloc1(nrecvs+1,&recv_waits);CHKERRQ(ierr); 108189d82c54SBarry Smith for (i=0; i<nrecvs; i++) { 108232dcc486SBarry Smith ierr = MPI_Irecv(recvs+2*nmax*i,2*nmax,MPIU_INT,MPI_ANY_SOURCE,tag1,comm,recv_waits+i);CHKERRQ(ierr); 108389d82c54SBarry Smith } 108489d82c54SBarry Smith 108589d82c54SBarry Smith /* pack messages containing lists of local nodes to owners */ 1086854ce69bSBarry Smith ierr = PetscMalloc1(2*n+1,&sends);CHKERRQ(ierr); 1087854ce69bSBarry Smith ierr = PetscMalloc1(size+1,&starts);CHKERRQ(ierr); 108889d82c54SBarry Smith starts[0] = 0; 1089f6e5521dSKarl Rupp for (i=1; i<size; i++) starts[i] = starts[i-1] + 2*nprocs[2*i-2]; 109089d82c54SBarry Smith for (i=0; i<n; i++) { 109189d82c54SBarry Smith sends[starts[owner[i]]++] = lindices[i]; 109230dcb7c9SBarry Smith sends[starts[owner[i]]++] = i; 109389d82c54SBarry Smith } 109489d82c54SBarry Smith ierr = PetscFree(owner);CHKERRQ(ierr); 109589d82c54SBarry Smith starts[0] = 0; 1096f6e5521dSKarl Rupp for (i=1; i<size; i++) starts[i] = starts[i-1] + 2*nprocs[2*i-2]; 109789d82c54SBarry Smith 109889d82c54SBarry Smith /* send the messages */ 1099854ce69bSBarry Smith ierr = PetscMalloc1(nsends+1,&send_waits);CHKERRQ(ierr); 1100854ce69bSBarry Smith ierr = PetscMalloc1(nsends+1,&dest);CHKERRQ(ierr); 110189d82c54SBarry Smith cnt = 0; 110289d82c54SBarry Smith for (i=0; i<size; i++) { 110327c402fcSBarry Smith if (nprocs[2*i]) { 110432dcc486SBarry Smith ierr = MPI_Isend(sends+starts[i],2*nprocs[2*i],MPIU_INT,i,tag1,comm,send_waits+cnt);CHKERRQ(ierr); 110530dcb7c9SBarry Smith dest[cnt] = i; 110689d82c54SBarry Smith cnt++; 110789d82c54SBarry Smith } 110889d82c54SBarry Smith } 110989d82c54SBarry Smith ierr = PetscFree(starts);CHKERRQ(ierr); 111089d82c54SBarry Smith 111189d82c54SBarry Smith /* wait on receives */ 1112854ce69bSBarry Smith ierr = PetscMalloc1(nrecvs+1,&source);CHKERRQ(ierr); 1113854ce69bSBarry Smith ierr = PetscMalloc1(nrecvs+1,&len);CHKERRQ(ierr); 111489d82c54SBarry Smith cnt = nrecvs; 1115854ce69bSBarry Smith ierr = PetscMalloc1(ng+1,&nownedsenders);CHKERRQ(ierr); 111632dcc486SBarry Smith ierr = PetscMemzero(nownedsenders,ng*sizeof(PetscInt));CHKERRQ(ierr); 111789d82c54SBarry Smith while (cnt) { 111889d82c54SBarry Smith ierr = MPI_Waitany(nrecvs,recv_waits,&imdex,&recv_status);CHKERRQ(ierr); 111989d82c54SBarry Smith /* unpack receives into our local space */ 112032dcc486SBarry Smith ierr = MPI_Get_count(&recv_status,MPIU_INT,&len[imdex]);CHKERRQ(ierr); 112189d82c54SBarry Smith source[imdex] = recv_status.MPI_SOURCE; 112230dcb7c9SBarry Smith len[imdex] = len[imdex]/2; 1123caba0dd0SBarry Smith /* count how many local owners for each of my global owned indices */ 112430dcb7c9SBarry Smith for (i=0; i<len[imdex]; i++) nownedsenders[recvs[2*imdex*nmax+2*i]-rstart]++; 112589d82c54SBarry Smith cnt--; 112689d82c54SBarry Smith } 112789d82c54SBarry Smith ierr = PetscFree(recv_waits);CHKERRQ(ierr); 112889d82c54SBarry Smith 112930dcb7c9SBarry Smith /* count how many globally owned indices are on an edge multiplied by how many processors own them. */ 1130bc8ff85bSBarry Smith nowned = 0; 1131bc8ff85bSBarry Smith nownedm = 0; 1132bc8ff85bSBarry Smith for (i=0; i<ng; i++) { 1133bc8ff85bSBarry Smith if (nownedsenders[i] > 1) {nownedm += nownedsenders[i]; nowned++;} 1134bc8ff85bSBarry Smith } 1135bc8ff85bSBarry Smith 1136bc8ff85bSBarry Smith /* create single array to contain rank of all local owners of each globally owned index */ 1137854ce69bSBarry Smith ierr = PetscMalloc1(nownedm+1,&ownedsenders);CHKERRQ(ierr); 1138854ce69bSBarry Smith ierr = PetscMalloc1(ng+1,&starts);CHKERRQ(ierr); 1139bc8ff85bSBarry Smith starts[0] = 0; 1140bc8ff85bSBarry Smith for (i=1; i<ng; i++) { 1141bc8ff85bSBarry Smith if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1]; 1142bc8ff85bSBarry Smith else starts[i] = starts[i-1]; 1143bc8ff85bSBarry Smith } 1144bc8ff85bSBarry Smith 114530dcb7c9SBarry Smith /* for each nontrival globally owned node list all arriving processors */ 1146bc8ff85bSBarry Smith for (i=0; i<nrecvs; i++) { 1147bc8ff85bSBarry Smith for (j=0; j<len[i]; j++) { 114830dcb7c9SBarry Smith node = recvs[2*i*nmax+2*j]-rstart; 1149f6e5521dSKarl Rupp if (nownedsenders[node] > 1) ownedsenders[starts[node]++] = source[i]; 1150bc8ff85bSBarry Smith } 1151bc8ff85bSBarry Smith } 1152bc8ff85bSBarry Smith 115307b52d57SBarry Smith if (debug) { /* ----------------------------------- */ 115430dcb7c9SBarry Smith starts[0] = 0; 115530dcb7c9SBarry Smith for (i=1; i<ng; i++) { 115630dcb7c9SBarry Smith if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1]; 115730dcb7c9SBarry Smith else starts[i] = starts[i-1]; 115830dcb7c9SBarry Smith } 115930dcb7c9SBarry Smith for (i=0; i<ng; i++) { 116030dcb7c9SBarry Smith if (nownedsenders[i] > 1) { 11617904a332SBarry Smith ierr = PetscSynchronizedPrintf(comm,"[%d] global node %D local owner processors: ",rank,i+rstart);CHKERRQ(ierr); 116230dcb7c9SBarry Smith for (j=0; j<nownedsenders[i]; j++) { 11637904a332SBarry Smith ierr = PetscSynchronizedPrintf(comm,"%D ",ownedsenders[starts[i]+j]);CHKERRQ(ierr); 116430dcb7c9SBarry Smith } 116530dcb7c9SBarry Smith ierr = PetscSynchronizedPrintf(comm,"\n");CHKERRQ(ierr); 116630dcb7c9SBarry Smith } 116730dcb7c9SBarry Smith } 11680ec8b6e3SBarry Smith ierr = PetscSynchronizedFlush(comm,PETSC_STDOUT);CHKERRQ(ierr); 116907b52d57SBarry Smith } /* ----------------------------------- */ 117030dcb7c9SBarry Smith 11713677ff5aSBarry Smith /* wait on original sends */ 11723a96401aSBarry Smith if (nsends) { 1173785e854fSJed Brown ierr = PetscMalloc1(nsends,&send_status);CHKERRQ(ierr); 11743a96401aSBarry Smith ierr = MPI_Waitall(nsends,send_waits,send_status);CHKERRQ(ierr); 11753a96401aSBarry Smith ierr = PetscFree(send_status);CHKERRQ(ierr); 11763a96401aSBarry Smith } 117789d82c54SBarry Smith ierr = PetscFree(send_waits);CHKERRQ(ierr); 11783a96401aSBarry Smith ierr = PetscFree(sends);CHKERRQ(ierr); 11793677ff5aSBarry Smith ierr = PetscFree(nprocs);CHKERRQ(ierr); 11803677ff5aSBarry Smith 11813677ff5aSBarry Smith /* pack messages to send back to local owners */ 118230dcb7c9SBarry Smith starts[0] = 0; 118330dcb7c9SBarry Smith for (i=1; i<ng; i++) { 118430dcb7c9SBarry Smith if (nownedsenders[i-1] > 1) starts[i] = starts[i-1] + nownedsenders[i-1]; 118530dcb7c9SBarry Smith else starts[i] = starts[i-1]; 118630dcb7c9SBarry Smith } 118730dcb7c9SBarry Smith nsends2 = nrecvs; 1188854ce69bSBarry Smith ierr = PetscMalloc1(nsends2+1,&nprocs);CHKERRQ(ierr); /* length of each message */ 118930dcb7c9SBarry Smith for (i=0; i<nrecvs; i++) { 119030dcb7c9SBarry Smith nprocs[i] = 1; 119130dcb7c9SBarry Smith for (j=0; j<len[i]; j++) { 119230dcb7c9SBarry Smith node = recvs[2*i*nmax+2*j]-rstart; 1193f6e5521dSKarl Rupp if (nownedsenders[node] > 1) nprocs[i] += 2 + nownedsenders[node]; 119430dcb7c9SBarry Smith } 119530dcb7c9SBarry Smith } 1196f6e5521dSKarl Rupp nt = 0; 1197f6e5521dSKarl Rupp for (i=0; i<nsends2; i++) nt += nprocs[i]; 1198f6e5521dSKarl Rupp 1199854ce69bSBarry Smith ierr = PetscMalloc1(nt+1,&sends2);CHKERRQ(ierr); 1200854ce69bSBarry Smith ierr = PetscMalloc1(nsends2+1,&starts2);CHKERRQ(ierr); 1201f6e5521dSKarl Rupp 1202f6e5521dSKarl Rupp starts2[0] = 0; 1203f6e5521dSKarl Rupp for (i=1; i<nsends2; i++) starts2[i] = starts2[i-1] + nprocs[i-1]; 120430dcb7c9SBarry Smith /* 120530dcb7c9SBarry Smith Each message is 1 + nprocs[i] long, and consists of 120630dcb7c9SBarry Smith (0) the number of nodes being sent back 120730dcb7c9SBarry Smith (1) the local node number, 120830dcb7c9SBarry Smith (2) the number of processors sharing it, 120930dcb7c9SBarry Smith (3) the processors sharing it 121030dcb7c9SBarry Smith */ 121130dcb7c9SBarry Smith for (i=0; i<nsends2; i++) { 121230dcb7c9SBarry Smith cnt = 1; 121330dcb7c9SBarry Smith sends2[starts2[i]] = 0; 121430dcb7c9SBarry Smith for (j=0; j<len[i]; j++) { 121530dcb7c9SBarry Smith node = recvs[2*i*nmax+2*j]-rstart; 121630dcb7c9SBarry Smith if (nownedsenders[node] > 1) { 121730dcb7c9SBarry Smith sends2[starts2[i]]++; 121830dcb7c9SBarry Smith sends2[starts2[i]+cnt++] = recvs[2*i*nmax+2*j+1]; 121930dcb7c9SBarry Smith sends2[starts2[i]+cnt++] = nownedsenders[node]; 122032dcc486SBarry Smith ierr = PetscMemcpy(&sends2[starts2[i]+cnt],&ownedsenders[starts[node]],nownedsenders[node]*sizeof(PetscInt));CHKERRQ(ierr); 122130dcb7c9SBarry Smith cnt += nownedsenders[node]; 122230dcb7c9SBarry Smith } 122330dcb7c9SBarry Smith } 122430dcb7c9SBarry Smith } 122530dcb7c9SBarry Smith 122630dcb7c9SBarry Smith /* receive the message lengths */ 122730dcb7c9SBarry Smith nrecvs2 = nsends; 1228854ce69bSBarry Smith ierr = PetscMalloc1(nrecvs2+1,&lens2);CHKERRQ(ierr); 1229854ce69bSBarry Smith ierr = PetscMalloc1(nrecvs2+1,&starts3);CHKERRQ(ierr); 1230854ce69bSBarry Smith ierr = PetscMalloc1(nrecvs2+1,&recv_waits);CHKERRQ(ierr); 123130dcb7c9SBarry Smith for (i=0; i<nrecvs2; i++) { 1232d44834fbSBarry Smith ierr = MPI_Irecv(&lens2[i],1,MPIU_INT,dest[i],tag2,comm,recv_waits+i);CHKERRQ(ierr); 123330dcb7c9SBarry Smith } 1234d44834fbSBarry Smith 12358a8e0b3aSBarry Smith /* send the message lengths */ 12368a8e0b3aSBarry Smith for (i=0; i<nsends2; i++) { 12378a8e0b3aSBarry Smith ierr = MPI_Send(&nprocs[i],1,MPIU_INT,source[i],tag2,comm);CHKERRQ(ierr); 12388a8e0b3aSBarry Smith } 12398a8e0b3aSBarry Smith 1240d44834fbSBarry Smith /* wait on receives of lens */ 12410c468ba9SBarry Smith if (nrecvs2) { 1242785e854fSJed Brown ierr = PetscMalloc1(nrecvs2,&recv_statuses);CHKERRQ(ierr); 1243d44834fbSBarry Smith ierr = MPI_Waitall(nrecvs2,recv_waits,recv_statuses);CHKERRQ(ierr); 1244d44834fbSBarry Smith ierr = PetscFree(recv_statuses);CHKERRQ(ierr); 12450c468ba9SBarry Smith } 1246a2ea699eSBarry Smith ierr = PetscFree(recv_waits);CHKERRQ(ierr); 1247d44834fbSBarry Smith 124830dcb7c9SBarry Smith starts3[0] = 0; 1249d44834fbSBarry Smith nt = 0; 125030dcb7c9SBarry Smith for (i=0; i<nrecvs2-1; i++) { 125130dcb7c9SBarry Smith starts3[i+1] = starts3[i] + lens2[i]; 1252d44834fbSBarry Smith nt += lens2[i]; 125330dcb7c9SBarry Smith } 125476466f69SStefano Zampini if (nrecvs2) nt += lens2[nrecvs2-1]; 1255d44834fbSBarry Smith 1256854ce69bSBarry Smith ierr = PetscMalloc1(nt+1,&recvs2);CHKERRQ(ierr); 1257854ce69bSBarry Smith ierr = PetscMalloc1(nrecvs2+1,&recv_waits);CHKERRQ(ierr); 125852b72c4aSBarry Smith for (i=0; i<nrecvs2; i++) { 125932dcc486SBarry Smith ierr = MPI_Irecv(recvs2+starts3[i],lens2[i],MPIU_INT,dest[i],tag3,comm,recv_waits+i);CHKERRQ(ierr); 126030dcb7c9SBarry Smith } 126130dcb7c9SBarry Smith 126230dcb7c9SBarry Smith /* send the messages */ 1263854ce69bSBarry Smith ierr = PetscMalloc1(nsends2+1,&send_waits);CHKERRQ(ierr); 126430dcb7c9SBarry Smith for (i=0; i<nsends2; i++) { 126532dcc486SBarry Smith ierr = MPI_Isend(sends2+starts2[i],nprocs[i],MPIU_INT,source[i],tag3,comm,send_waits+i);CHKERRQ(ierr); 126630dcb7c9SBarry Smith } 126730dcb7c9SBarry Smith 126830dcb7c9SBarry Smith /* wait on receives */ 12690c468ba9SBarry Smith if (nrecvs2) { 1270785e854fSJed Brown ierr = PetscMalloc1(nrecvs2,&recv_statuses);CHKERRQ(ierr); 127130dcb7c9SBarry Smith ierr = MPI_Waitall(nrecvs2,recv_waits,recv_statuses);CHKERRQ(ierr); 127230dcb7c9SBarry Smith ierr = PetscFree(recv_statuses);CHKERRQ(ierr); 12730c468ba9SBarry Smith } 127430dcb7c9SBarry Smith ierr = PetscFree(recv_waits);CHKERRQ(ierr); 127530dcb7c9SBarry Smith ierr = PetscFree(nprocs);CHKERRQ(ierr); 127630dcb7c9SBarry Smith 127707b52d57SBarry Smith if (debug) { /* ----------------------------------- */ 127830dcb7c9SBarry Smith cnt = 0; 127930dcb7c9SBarry Smith for (i=0; i<nrecvs2; i++) { 128030dcb7c9SBarry Smith nt = recvs2[cnt++]; 128130dcb7c9SBarry Smith for (j=0; j<nt; j++) { 12827904a332SBarry Smith ierr = PetscSynchronizedPrintf(comm,"[%d] local node %D number of subdomains %D: ",rank,recvs2[cnt],recvs2[cnt+1]);CHKERRQ(ierr); 128330dcb7c9SBarry Smith for (k=0; k<recvs2[cnt+1]; k++) { 12847904a332SBarry Smith ierr = PetscSynchronizedPrintf(comm,"%D ",recvs2[cnt+2+k]);CHKERRQ(ierr); 128530dcb7c9SBarry Smith } 128630dcb7c9SBarry Smith cnt += 2 + recvs2[cnt+1]; 128730dcb7c9SBarry Smith ierr = PetscSynchronizedPrintf(comm,"\n");CHKERRQ(ierr); 128830dcb7c9SBarry Smith } 128930dcb7c9SBarry Smith } 12900ec8b6e3SBarry Smith ierr = PetscSynchronizedFlush(comm,PETSC_STDOUT);CHKERRQ(ierr); 129107b52d57SBarry Smith } /* ----------------------------------- */ 129230dcb7c9SBarry Smith 129330dcb7c9SBarry Smith /* count number subdomains for each local node */ 1294785e854fSJed Brown ierr = PetscMalloc1(size,&nprocs);CHKERRQ(ierr); 129532dcc486SBarry Smith ierr = PetscMemzero(nprocs,size*sizeof(PetscInt));CHKERRQ(ierr); 129630dcb7c9SBarry Smith cnt = 0; 129730dcb7c9SBarry Smith for (i=0; i<nrecvs2; i++) { 129830dcb7c9SBarry Smith nt = recvs2[cnt++]; 129930dcb7c9SBarry Smith for (j=0; j<nt; j++) { 1300f6e5521dSKarl Rupp for (k=0; k<recvs2[cnt+1]; k++) nprocs[recvs2[cnt+2+k]]++; 130130dcb7c9SBarry Smith cnt += 2 + recvs2[cnt+1]; 130230dcb7c9SBarry Smith } 130330dcb7c9SBarry Smith } 130430dcb7c9SBarry Smith nt = 0; for (i=0; i<size; i++) nt += (nprocs[i] > 0); 130530dcb7c9SBarry Smith *nproc = nt; 1306854ce69bSBarry Smith ierr = PetscMalloc1(nt+1,procs);CHKERRQ(ierr); 1307854ce69bSBarry Smith ierr = PetscMalloc1(nt+1,numprocs);CHKERRQ(ierr); 1308854ce69bSBarry Smith ierr = PetscMalloc1(nt+1,indices);CHKERRQ(ierr); 13090298fd71SBarry Smith for (i=0;i<nt+1;i++) (*indices)[i]=NULL; 1310785e854fSJed Brown ierr = PetscMalloc1(size,&bprocs);CHKERRQ(ierr); 131130dcb7c9SBarry Smith cnt = 0; 131230dcb7c9SBarry Smith for (i=0; i<size; i++) { 131330dcb7c9SBarry Smith if (nprocs[i] > 0) { 131430dcb7c9SBarry Smith bprocs[i] = cnt; 131530dcb7c9SBarry Smith (*procs)[cnt] = i; 131630dcb7c9SBarry Smith (*numprocs)[cnt] = nprocs[i]; 1317785e854fSJed Brown ierr = PetscMalloc1(nprocs[i],&(*indices)[cnt]);CHKERRQ(ierr); 131830dcb7c9SBarry Smith cnt++; 131930dcb7c9SBarry Smith } 132030dcb7c9SBarry Smith } 132130dcb7c9SBarry Smith 132230dcb7c9SBarry Smith /* make the list of subdomains for each nontrivial local node */ 132332dcc486SBarry Smith ierr = PetscMemzero(*numprocs,nt*sizeof(PetscInt));CHKERRQ(ierr); 132430dcb7c9SBarry Smith cnt = 0; 132530dcb7c9SBarry Smith for (i=0; i<nrecvs2; i++) { 132630dcb7c9SBarry Smith nt = recvs2[cnt++]; 132730dcb7c9SBarry Smith for (j=0; j<nt; j++) { 1328f6e5521dSKarl Rupp for (k=0; k<recvs2[cnt+1]; k++) (*indices)[bprocs[recvs2[cnt+2+k]]][(*numprocs)[bprocs[recvs2[cnt+2+k]]]++] = recvs2[cnt]; 132930dcb7c9SBarry Smith cnt += 2 + recvs2[cnt+1]; 133030dcb7c9SBarry Smith } 133130dcb7c9SBarry Smith } 133230dcb7c9SBarry Smith ierr = PetscFree(bprocs);CHKERRQ(ierr); 133307b52d57SBarry Smith ierr = PetscFree(recvs2);CHKERRQ(ierr); 133430dcb7c9SBarry Smith 133507b52d57SBarry Smith /* sort the node indexing by their global numbers */ 133607b52d57SBarry Smith nt = *nproc; 133707b52d57SBarry Smith for (i=0; i<nt; i++) { 1338854ce69bSBarry Smith ierr = PetscMalloc1((*numprocs)[i],&tmp);CHKERRQ(ierr); 1339f6e5521dSKarl Rupp for (j=0; j<(*numprocs)[i]; j++) tmp[j] = lindices[(*indices)[i][j]]; 134007b52d57SBarry Smith ierr = PetscSortIntWithArray((*numprocs)[i],tmp,(*indices)[i]);CHKERRQ(ierr); 134107b52d57SBarry Smith ierr = PetscFree(tmp);CHKERRQ(ierr); 134207b52d57SBarry Smith } 134307b52d57SBarry Smith 134407b52d57SBarry Smith if (debug) { /* ----------------------------------- */ 134530dcb7c9SBarry Smith nt = *nproc; 134630dcb7c9SBarry Smith for (i=0; i<nt; i++) { 13477904a332SBarry Smith ierr = PetscSynchronizedPrintf(comm,"[%d] subdomain %D number of indices %D: ",rank,(*procs)[i],(*numprocs)[i]);CHKERRQ(ierr); 134830dcb7c9SBarry Smith for (j=0; j<(*numprocs)[i]; j++) { 13497904a332SBarry Smith ierr = PetscSynchronizedPrintf(comm,"%D ",(*indices)[i][j]);CHKERRQ(ierr); 135030dcb7c9SBarry Smith } 135130dcb7c9SBarry Smith ierr = PetscSynchronizedPrintf(comm,"\n");CHKERRQ(ierr); 135230dcb7c9SBarry Smith } 13530ec8b6e3SBarry Smith ierr = PetscSynchronizedFlush(comm,PETSC_STDOUT);CHKERRQ(ierr); 135407b52d57SBarry Smith } /* ----------------------------------- */ 135530dcb7c9SBarry Smith 135630dcb7c9SBarry Smith /* wait on sends */ 135730dcb7c9SBarry Smith if (nsends2) { 1358785e854fSJed Brown ierr = PetscMalloc1(nsends2,&send_status);CHKERRQ(ierr); 135930dcb7c9SBarry Smith ierr = MPI_Waitall(nsends2,send_waits,send_status);CHKERRQ(ierr); 136030dcb7c9SBarry Smith ierr = PetscFree(send_status);CHKERRQ(ierr); 136130dcb7c9SBarry Smith } 136230dcb7c9SBarry Smith 136330dcb7c9SBarry Smith ierr = PetscFree(starts3);CHKERRQ(ierr); 136430dcb7c9SBarry Smith ierr = PetscFree(dest);CHKERRQ(ierr); 136530dcb7c9SBarry Smith ierr = PetscFree(send_waits);CHKERRQ(ierr); 13663677ff5aSBarry Smith 1367bc8ff85bSBarry Smith ierr = PetscFree(nownedsenders);CHKERRQ(ierr); 1368bc8ff85bSBarry Smith ierr = PetscFree(ownedsenders);CHKERRQ(ierr); 1369bc8ff85bSBarry Smith ierr = PetscFree(starts);CHKERRQ(ierr); 137030dcb7c9SBarry Smith ierr = PetscFree(starts2);CHKERRQ(ierr); 137130dcb7c9SBarry Smith ierr = PetscFree(lens2);CHKERRQ(ierr); 137289d82c54SBarry Smith 137389d82c54SBarry Smith ierr = PetscFree(source);CHKERRQ(ierr); 137497f1f81fSBarry Smith ierr = PetscFree(len);CHKERRQ(ierr); 137589d82c54SBarry Smith ierr = PetscFree(recvs);CHKERRQ(ierr); 13763a96401aSBarry Smith ierr = PetscFree(nprocs);CHKERRQ(ierr); 137730dcb7c9SBarry Smith ierr = PetscFree(sends2);CHKERRQ(ierr); 137824cf384cSBarry Smith 137924cf384cSBarry Smith /* put the information about myself as the first entry in the list */ 138024cf384cSBarry Smith first_procs = (*procs)[0]; 138124cf384cSBarry Smith first_numprocs = (*numprocs)[0]; 138224cf384cSBarry Smith first_indices = (*indices)[0]; 138324cf384cSBarry Smith for (i=0; i<*nproc; i++) { 138424cf384cSBarry Smith if ((*procs)[i] == rank) { 138524cf384cSBarry Smith (*procs)[0] = (*procs)[i]; 138624cf384cSBarry Smith (*numprocs)[0] = (*numprocs)[i]; 138724cf384cSBarry Smith (*indices)[0] = (*indices)[i]; 138824cf384cSBarry Smith (*procs)[i] = first_procs; 138924cf384cSBarry Smith (*numprocs)[i] = first_numprocs; 139024cf384cSBarry Smith (*indices)[i] = first_indices; 139124cf384cSBarry Smith break; 139224cf384cSBarry Smith } 139324cf384cSBarry Smith } 1394268a049cSStefano Zampini 1395268a049cSStefano Zampini /* save info for reuse */ 1396268a049cSStefano Zampini mapping->info_nproc = *nproc; 1397268a049cSStefano Zampini mapping->info_procs = *procs; 1398268a049cSStefano Zampini mapping->info_numprocs = *numprocs; 1399268a049cSStefano Zampini mapping->info_indices = *indices; 1400268a049cSStefano Zampini mapping->info_cached = PETSC_TRUE; 140189d82c54SBarry Smith PetscFunctionReturn(0); 140289d82c54SBarry Smith } 140389d82c54SBarry Smith 14046a818285SBarry Smith /*@C 14056a818285SBarry Smith ISLocalToGlobalMappingRestoreBlockInfo - Frees the memory allocated by ISLocalToGlobalMappingGetBlockInfo() 14066a818285SBarry Smith 14076a818285SBarry Smith Collective on ISLocalToGlobalMapping 14086a818285SBarry Smith 14096a818285SBarry Smith Input Parameters: 14106a818285SBarry Smith . mapping - the mapping from local to global indexing 14116a818285SBarry Smith 14126a818285SBarry Smith Output Parameter: 14136a818285SBarry Smith + nproc - number of processors that are connected to this one 14146a818285SBarry Smith . proc - neighboring processors 14156a818285SBarry Smith . numproc - number of indices for each processor 14166a818285SBarry Smith - indices - indices of local nodes shared with neighbor (sorted by global numbering) 14176a818285SBarry Smith 14186a818285SBarry Smith Level: advanced 14196a818285SBarry Smith 14206a818285SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(), 14216a818285SBarry Smith ISLocalToGlobalMappingGetInfo() 14226a818285SBarry Smith @*/ 14236a818285SBarry Smith PetscErrorCode ISLocalToGlobalMappingRestoreBlockInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[]) 14246a818285SBarry Smith { 14256a818285SBarry Smith PetscErrorCode ierr; 14266a818285SBarry Smith 14276a818285SBarry Smith PetscFunctionBegin; 1428cbc1caf0SMatthew G. Knepley PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1); 1429268a049cSStefano Zampini if (mapping->info_free) { 14306a818285SBarry Smith ierr = PetscFree(*numprocs);CHKERRQ(ierr); 14316a818285SBarry Smith if (*indices) { 1432268a049cSStefano Zampini PetscInt i; 1433268a049cSStefano Zampini 14346a818285SBarry Smith ierr = PetscFree((*indices)[0]);CHKERRQ(ierr); 14356a818285SBarry Smith for (i=1; i<*nproc; i++) { 14366a818285SBarry Smith ierr = PetscFree((*indices)[i]);CHKERRQ(ierr); 14376a818285SBarry Smith } 14386a818285SBarry Smith ierr = PetscFree(*indices);CHKERRQ(ierr); 14396a818285SBarry Smith } 1440268a049cSStefano Zampini } 1441268a049cSStefano Zampini *nproc = 0; 1442268a049cSStefano Zampini *procs = NULL; 1443268a049cSStefano Zampini *numprocs = NULL; 1444268a049cSStefano Zampini *indices = NULL; 14456a818285SBarry Smith PetscFunctionReturn(0); 14466a818285SBarry Smith } 14476a818285SBarry Smith 14486a818285SBarry Smith /*@C 14496a818285SBarry Smith ISLocalToGlobalMappingGetInfo - Gets the neighbor information for each processor and 14506a818285SBarry Smith each index shared by more than one processor 14516a818285SBarry Smith 14526a818285SBarry Smith Collective on ISLocalToGlobalMapping 14536a818285SBarry Smith 14546a818285SBarry Smith Input Parameters: 14556a818285SBarry Smith . mapping - the mapping from local to global indexing 14566a818285SBarry Smith 14576a818285SBarry Smith Output Parameter: 14586a818285SBarry Smith + nproc - number of processors that are connected to this one 14596a818285SBarry Smith . proc - neighboring processors 14606a818285SBarry Smith . numproc - number of indices for each subdomain (processor) 14616a818285SBarry Smith - indices - indices of nodes (in local numbering) shared with neighbors (sorted by global numbering) 14626a818285SBarry Smith 14636a818285SBarry Smith Level: advanced 14646a818285SBarry Smith 14656a818285SBarry Smith Concepts: mapping^local to global 14666a818285SBarry Smith 14671bd0b88eSStefano Zampini Notes: The user needs to call ISLocalToGlobalMappingRestoreInfo when the data is no longer needed. 14681bd0b88eSStefano Zampini 14696a818285SBarry Smith Fortran Usage: 14706a818285SBarry Smith $ ISLocalToGlobalMpngGetInfoSize(ISLocalToGlobalMapping,PetscInt nproc,PetscInt numprocmax,ierr) followed by 14716a818285SBarry Smith $ ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping,PetscInt nproc, PetscInt procs[nproc],PetscInt numprocs[nproc], 14726a818285SBarry Smith PetscInt indices[nproc][numprocmax],ierr) 14736a818285SBarry Smith There is no ISLocalToGlobalMappingRestoreInfo() in Fortran. You must make sure that procs[], numprocs[] and 14746a818285SBarry Smith indices[][] are large enough arrays, either by allocating them dynamically or defining static ones large enough. 14756a818285SBarry Smith 14766a818285SBarry Smith 14776a818285SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(), 14786a818285SBarry Smith ISLocalToGlobalMappingRestoreInfo() 14796a818285SBarry Smith @*/ 14806a818285SBarry Smith PetscErrorCode ISLocalToGlobalMappingGetInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[]) 14816a818285SBarry Smith { 14826a818285SBarry Smith PetscErrorCode ierr; 1483268a049cSStefano Zampini PetscInt **bindices = NULL,*bnumprocs = NULL,bs = mapping->bs,i,j,k; 14846a818285SBarry Smith 14856a818285SBarry Smith PetscFunctionBegin; 14866a818285SBarry Smith PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1); 1487268a049cSStefano Zampini ierr = ISLocalToGlobalMappingGetBlockInfo(mapping,nproc,procs,&bnumprocs,&bindices);CHKERRQ(ierr); 1488268a049cSStefano Zampini if (bs > 1) { /* we need to expand the cached info */ 1489732f65e3SBarry Smith ierr = PetscCalloc1(*nproc,&*indices);CHKERRQ(ierr); 1490268a049cSStefano Zampini ierr = PetscCalloc1(*nproc,&*numprocs);CHKERRQ(ierr); 14916a818285SBarry Smith for (i=0; i<*nproc; i++) { 1492268a049cSStefano Zampini ierr = PetscMalloc1(bs*bnumprocs[i],&(*indices)[i]);CHKERRQ(ierr); 1493268a049cSStefano Zampini for (j=0; j<bnumprocs[i]; j++) { 14946a818285SBarry Smith for (k=0; k<bs; k++) { 14956a818285SBarry Smith (*indices)[i][j*bs+k] = bs*bindices[i][j] + k; 14966a818285SBarry Smith } 14976a818285SBarry Smith } 1498268a049cSStefano Zampini (*numprocs)[i] = bnumprocs[i]*bs; 14996a818285SBarry Smith } 1500268a049cSStefano Zampini mapping->info_free = PETSC_TRUE; 1501268a049cSStefano Zampini } else { 1502268a049cSStefano Zampini *numprocs = bnumprocs; 1503268a049cSStefano Zampini *indices = bindices; 15046a818285SBarry Smith } 15056a818285SBarry Smith PetscFunctionReturn(0); 15066a818285SBarry Smith } 15076a818285SBarry Smith 150807b52d57SBarry Smith /*@C 150907b52d57SBarry Smith ISLocalToGlobalMappingRestoreInfo - Frees the memory allocated by ISLocalToGlobalMappingGetInfo() 151089d82c54SBarry Smith 151107b52d57SBarry Smith Collective on ISLocalToGlobalMapping 151207b52d57SBarry Smith 151307b52d57SBarry Smith Input Parameters: 151407b52d57SBarry Smith . mapping - the mapping from local to global indexing 151507b52d57SBarry Smith 151607b52d57SBarry Smith Output Parameter: 151707b52d57SBarry Smith + nproc - number of processors that are connected to this one 151807b52d57SBarry Smith . proc - neighboring processors 151907b52d57SBarry Smith . numproc - number of indices for each processor 152007b52d57SBarry Smith - indices - indices of local nodes shared with neighbor (sorted by global numbering) 152107b52d57SBarry Smith 152207b52d57SBarry Smith Level: advanced 152307b52d57SBarry Smith 152407b52d57SBarry Smith .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(), 152507b52d57SBarry Smith ISLocalToGlobalMappingGetInfo() 152607b52d57SBarry Smith @*/ 15277087cfbeSBarry Smith PetscErrorCode ISLocalToGlobalMappingRestoreInfo(ISLocalToGlobalMapping mapping,PetscInt *nproc,PetscInt *procs[],PetscInt *numprocs[],PetscInt **indices[]) 152807b52d57SBarry Smith { 15296849ba73SBarry Smith PetscErrorCode ierr; 153007b52d57SBarry Smith 153107b52d57SBarry Smith PetscFunctionBegin; 15326a818285SBarry Smith ierr = ISLocalToGlobalMappingRestoreBlockInfo(mapping,nproc,procs,numprocs,indices);CHKERRQ(ierr); 153307b52d57SBarry Smith PetscFunctionReturn(0); 153407b52d57SBarry Smith } 153586994e45SJed Brown 153686994e45SJed Brown /*@C 15371bd0b88eSStefano Zampini ISLocalToGlobalMappingGetNodeInfo - Gets the neighbor information for each node 15381bd0b88eSStefano Zampini 15391bd0b88eSStefano Zampini Collective on ISLocalToGlobalMapping 15401bd0b88eSStefano Zampini 15411bd0b88eSStefano Zampini Input Parameters: 15421bd0b88eSStefano Zampini . mapping - the mapping from local to global indexing 15431bd0b88eSStefano Zampini 15441bd0b88eSStefano Zampini Output Parameter: 15451bd0b88eSStefano Zampini + nnodes - number of local nodes (same ISLocalToGlobalMappingGetSize()) 15461bd0b88eSStefano Zampini . count - number of neighboring processors per node 15471bd0b88eSStefano Zampini - indices - indices of processes sharing the node (sorted) 15481bd0b88eSStefano Zampini 15491bd0b88eSStefano Zampini Level: advanced 15501bd0b88eSStefano Zampini 15511bd0b88eSStefano Zampini Concepts: mapping^local to global 15521bd0b88eSStefano Zampini 15531bd0b88eSStefano Zampini Notes: The user needs to call ISLocalToGlobalMappingRestoreInfo when the data is no longer needed. 15541bd0b88eSStefano Zampini 15551bd0b88eSStefano Zampini .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(), 15561bd0b88eSStefano Zampini ISLocalToGlobalMappingGetInfo(), ISLocalToGlobalMappingRestoreNodeInfo() 15571bd0b88eSStefano Zampini @*/ 15581bd0b88eSStefano Zampini PetscErrorCode ISLocalToGlobalMappingGetNodeInfo(ISLocalToGlobalMapping mapping,PetscInt *nnodes,PetscInt *count[],PetscInt **indices[]) 15591bd0b88eSStefano Zampini { 15601bd0b88eSStefano Zampini PetscInt n; 15611bd0b88eSStefano Zampini PetscErrorCode ierr; 15621bd0b88eSStefano Zampini 15631bd0b88eSStefano Zampini PetscFunctionBegin; 15641bd0b88eSStefano Zampini PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1); 15651bd0b88eSStefano Zampini ierr = ISLocalToGlobalMappingGetSize(mapping,&n);CHKERRQ(ierr); 15661bd0b88eSStefano Zampini if (!mapping->info_nodec) { 15671bd0b88eSStefano Zampini PetscInt i,m,n_neigh,*neigh,*n_shared,**shared; 15681bd0b88eSStefano Zampini 15691bd0b88eSStefano Zampini ierr = PetscCalloc1(n+1,&mapping->info_nodec);CHKERRQ(ierr); /* always allocate to flag setup */ 15701bd0b88eSStefano Zampini ierr = PetscMalloc1(n,&mapping->info_nodei);CHKERRQ(ierr); 15711bd0b88eSStefano Zampini ierr = ISLocalToGlobalMappingGetInfo(mapping,&n_neigh,&neigh,&n_shared,&shared);CHKERRQ(ierr); 15721bd0b88eSStefano Zampini for (i=0,m=0;i<n;i++) { mapping->info_nodec[i] = 1; m++; } 15731bd0b88eSStefano Zampini for (i=1;i<n_neigh;i++) { 15741bd0b88eSStefano Zampini PetscInt j; 15751bd0b88eSStefano Zampini 15761bd0b88eSStefano Zampini m += n_shared[i]; 15771bd0b88eSStefano Zampini for (j=0;j<n_shared[i];j++) mapping->info_nodec[shared[i][j]] += 1; 15781bd0b88eSStefano Zampini } 15791bd0b88eSStefano Zampini if (n) { ierr = PetscMalloc1(m,&mapping->info_nodei[0]);CHKERRQ(ierr); } 15801bd0b88eSStefano Zampini for (i=1;i<n;i++) mapping->info_nodei[i] = mapping->info_nodei[i-1] + mapping->info_nodec[i-1]; 15811bd0b88eSStefano Zampini ierr = PetscMemzero(mapping->info_nodec,n*sizeof(PetscInt));CHKERRQ(ierr); 15821bd0b88eSStefano Zampini for (i=0;i<n;i++) { mapping->info_nodec[i] = 1; mapping->info_nodei[i][0] = neigh[0]; } 15831bd0b88eSStefano Zampini for (i=1;i<n_neigh;i++) { 15841bd0b88eSStefano Zampini PetscInt j; 15851bd0b88eSStefano Zampini 15861bd0b88eSStefano Zampini for (j=0;j<n_shared[i];j++) { 15871bd0b88eSStefano Zampini PetscInt k = shared[i][j]; 15881bd0b88eSStefano Zampini 15891bd0b88eSStefano Zampini mapping->info_nodei[k][mapping->info_nodec[k]] = neigh[i]; 15901bd0b88eSStefano Zampini mapping->info_nodec[k] += 1; 15911bd0b88eSStefano Zampini } 15921bd0b88eSStefano Zampini } 15931bd0b88eSStefano Zampini for (i=0;i<n;i++) { ierr = PetscSortRemoveDupsInt(&mapping->info_nodec[i],mapping->info_nodei[i]);CHKERRQ(ierr); } 15941bd0b88eSStefano Zampini ierr = ISLocalToGlobalMappingRestoreInfo(mapping,&n_neigh,&neigh,&n_shared,&shared);CHKERRQ(ierr); 15951bd0b88eSStefano Zampini } 15961bd0b88eSStefano Zampini if (nnodes) *nnodes = n; 15971bd0b88eSStefano Zampini if (count) *count = mapping->info_nodec; 15981bd0b88eSStefano Zampini if (indices) *indices = mapping->info_nodei; 15991bd0b88eSStefano Zampini PetscFunctionReturn(0); 16001bd0b88eSStefano Zampini } 16011bd0b88eSStefano Zampini 16021bd0b88eSStefano Zampini /*@C 16031bd0b88eSStefano Zampini ISLocalToGlobalMappingRestoreNodeInfo - Frees the memory allocated by ISLocalToGlobalMappingGetNodeInfo() 16041bd0b88eSStefano Zampini 16051bd0b88eSStefano Zampini Collective on ISLocalToGlobalMapping 16061bd0b88eSStefano Zampini 16071bd0b88eSStefano Zampini Input Parameters: 16081bd0b88eSStefano Zampini . mapping - the mapping from local to global indexing 16091bd0b88eSStefano Zampini 16101bd0b88eSStefano Zampini Output Parameter: 16111bd0b88eSStefano Zampini + nnodes - number of local nodes 16121bd0b88eSStefano Zampini . count - number of neighboring processors per node 16131bd0b88eSStefano Zampini - indices - indices of processes sharing the node (sorted) 16141bd0b88eSStefano Zampini 16151bd0b88eSStefano Zampini Level: advanced 16161bd0b88eSStefano Zampini 16171bd0b88eSStefano Zampini .seealso: ISLocalToGlobalMappingDestroy(), ISLocalToGlobalMappingCreateIS(), ISLocalToGlobalMappingCreate(), 16181bd0b88eSStefano Zampini ISLocalToGlobalMappingGetInfo() 16191bd0b88eSStefano Zampini @*/ 16201bd0b88eSStefano Zampini PetscErrorCode ISLocalToGlobalMappingRestoreNodeInfo(ISLocalToGlobalMapping mapping,PetscInt *nnodes,PetscInt *count[],PetscInt **indices[]) 16211bd0b88eSStefano Zampini { 16221bd0b88eSStefano Zampini PetscFunctionBegin; 16231bd0b88eSStefano Zampini PetscValidHeaderSpecific(mapping,IS_LTOGM_CLASSID,1); 16241bd0b88eSStefano Zampini if (nnodes) *nnodes = 0; 16251bd0b88eSStefano Zampini if (count) *count = NULL; 16261bd0b88eSStefano Zampini if (indices) *indices = NULL; 16271bd0b88eSStefano Zampini PetscFunctionReturn(0); 16281bd0b88eSStefano Zampini } 16291bd0b88eSStefano Zampini 16301bd0b88eSStefano Zampini /*@C 1631107e9a97SBarry Smith ISLocalToGlobalMappingGetIndices - Get global indices for every local point that is mapped 163286994e45SJed Brown 163386994e45SJed Brown Not Collective 163486994e45SJed Brown 163586994e45SJed Brown Input Arguments: 163686994e45SJed Brown . ltog - local to global mapping 163786994e45SJed Brown 163886994e45SJed Brown Output Arguments: 1639565245c5SBarry Smith . array - array of indices, the length of this array may be obtained with ISLocalToGlobalMappingGetSize() 164086994e45SJed Brown 164186994e45SJed Brown Level: advanced 164286994e45SJed Brown 164395452b02SPatrick Sanan Notes: 164495452b02SPatrick Sanan ISLocalToGlobalMappingGetSize() returns the length the this array 1645107e9a97SBarry Smith 1646107e9a97SBarry Smith .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingRestoreIndices(), ISLocalToGlobalMappingGetBlockIndices(), ISLocalToGlobalMappingRestoreBlockIndices() 164786994e45SJed Brown @*/ 16487087cfbeSBarry Smith PetscErrorCode ISLocalToGlobalMappingGetIndices(ISLocalToGlobalMapping ltog,const PetscInt **array) 164986994e45SJed Brown { 165086994e45SJed Brown PetscFunctionBegin; 165186994e45SJed Brown PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1); 165286994e45SJed Brown PetscValidPointer(array,2); 165345b6f7e9SBarry Smith if (ltog->bs == 1) { 165486994e45SJed Brown *array = ltog->indices; 165545b6f7e9SBarry Smith } else { 165645b6f7e9SBarry Smith PetscInt *jj,k,i,j,n = ltog->n, bs = ltog->bs; 165745b6f7e9SBarry Smith const PetscInt *ii; 165845b6f7e9SBarry Smith PetscErrorCode ierr; 165945b6f7e9SBarry Smith 166045b6f7e9SBarry Smith ierr = PetscMalloc1(bs*n,&jj);CHKERRQ(ierr); 166145b6f7e9SBarry Smith *array = jj; 166245b6f7e9SBarry Smith k = 0; 166345b6f7e9SBarry Smith ii = ltog->indices; 166445b6f7e9SBarry Smith for (i=0; i<n; i++) 166545b6f7e9SBarry Smith for (j=0; j<bs; j++) 166645b6f7e9SBarry Smith jj[k++] = bs*ii[i] + j; 166745b6f7e9SBarry Smith } 166886994e45SJed Brown PetscFunctionReturn(0); 166986994e45SJed Brown } 167086994e45SJed Brown 167186994e45SJed Brown /*@C 1672193a2b41SJulian Andrej ISLocalToGlobalMappingRestoreIndices - Restore indices obtained with ISLocalToGlobalMappingGetIndices() 167386994e45SJed Brown 167486994e45SJed Brown Not Collective 167586994e45SJed Brown 167686994e45SJed Brown Input Arguments: 167786994e45SJed Brown + ltog - local to global mapping 167886994e45SJed Brown - array - array of indices 167986994e45SJed Brown 168086994e45SJed Brown Level: advanced 168186994e45SJed Brown 168286994e45SJed Brown .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingGetIndices() 168386994e45SJed Brown @*/ 16847087cfbeSBarry Smith PetscErrorCode ISLocalToGlobalMappingRestoreIndices(ISLocalToGlobalMapping ltog,const PetscInt **array) 168586994e45SJed Brown { 168686994e45SJed Brown PetscFunctionBegin; 168786994e45SJed Brown PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1); 168886994e45SJed Brown PetscValidPointer(array,2); 168945b6f7e9SBarry Smith if (ltog->bs == 1 && *array != ltog->indices) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Trying to return mismatched pointer"); 169045b6f7e9SBarry Smith 169145b6f7e9SBarry Smith if (ltog->bs > 1) { 169245b6f7e9SBarry Smith PetscErrorCode ierr; 169345b6f7e9SBarry Smith ierr = PetscFree(*(void**)array);CHKERRQ(ierr); 169445b6f7e9SBarry Smith } 169545b6f7e9SBarry Smith PetscFunctionReturn(0); 169645b6f7e9SBarry Smith } 169745b6f7e9SBarry Smith 169845b6f7e9SBarry Smith /*@C 169945b6f7e9SBarry Smith ISLocalToGlobalMappingGetBlockIndices - Get global indices for every local block 170045b6f7e9SBarry Smith 170145b6f7e9SBarry Smith Not Collective 170245b6f7e9SBarry Smith 170345b6f7e9SBarry Smith Input Arguments: 170445b6f7e9SBarry Smith . ltog - local to global mapping 170545b6f7e9SBarry Smith 170645b6f7e9SBarry Smith Output Arguments: 170745b6f7e9SBarry Smith . array - array of indices 170845b6f7e9SBarry Smith 170945b6f7e9SBarry Smith Level: advanced 171045b6f7e9SBarry Smith 171145b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingRestoreBlockIndices() 171245b6f7e9SBarry Smith @*/ 171345b6f7e9SBarry Smith PetscErrorCode ISLocalToGlobalMappingGetBlockIndices(ISLocalToGlobalMapping ltog,const PetscInt **array) 171445b6f7e9SBarry Smith { 171545b6f7e9SBarry Smith PetscFunctionBegin; 171645b6f7e9SBarry Smith PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1); 171745b6f7e9SBarry Smith PetscValidPointer(array,2); 171845b6f7e9SBarry Smith *array = ltog->indices; 171945b6f7e9SBarry Smith PetscFunctionReturn(0); 172045b6f7e9SBarry Smith } 172145b6f7e9SBarry Smith 172245b6f7e9SBarry Smith /*@C 172345b6f7e9SBarry Smith ISLocalToGlobalMappingRestoreBlockIndices - Restore indices obtained with ISLocalToGlobalMappingGetBlockIndices() 172445b6f7e9SBarry Smith 172545b6f7e9SBarry Smith Not Collective 172645b6f7e9SBarry Smith 172745b6f7e9SBarry Smith Input Arguments: 172845b6f7e9SBarry Smith + ltog - local to global mapping 172945b6f7e9SBarry Smith - array - array of indices 173045b6f7e9SBarry Smith 173145b6f7e9SBarry Smith Level: advanced 173245b6f7e9SBarry Smith 173345b6f7e9SBarry Smith .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingApply(), ISLocalToGlobalMappingGetIndices() 173445b6f7e9SBarry Smith @*/ 173545b6f7e9SBarry Smith PetscErrorCode ISLocalToGlobalMappingRestoreBlockIndices(ISLocalToGlobalMapping ltog,const PetscInt **array) 173645b6f7e9SBarry Smith { 173745b6f7e9SBarry Smith PetscFunctionBegin; 173845b6f7e9SBarry Smith PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1); 173945b6f7e9SBarry Smith PetscValidPointer(array,2); 174086994e45SJed Brown if (*array != ltog->indices) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_BADPTR,"Trying to return mismatched pointer"); 17410298fd71SBarry Smith *array = NULL; 174286994e45SJed Brown PetscFunctionReturn(0); 174386994e45SJed Brown } 1744f7efa3c7SJed Brown 1745f7efa3c7SJed Brown /*@C 1746f7efa3c7SJed Brown ISLocalToGlobalMappingConcatenate - Create a new mapping that concatenates a list of mappings 1747f7efa3c7SJed Brown 1748f7efa3c7SJed Brown Not Collective 1749f7efa3c7SJed Brown 1750f7efa3c7SJed Brown Input Arguments: 1751f7efa3c7SJed Brown + comm - communicator for the new mapping, must contain the communicator of every mapping to concatenate 1752f7efa3c7SJed Brown . n - number of mappings to concatenate 1753f7efa3c7SJed Brown - ltogs - local to global mappings 1754f7efa3c7SJed Brown 1755f7efa3c7SJed Brown Output Arguments: 1756f7efa3c7SJed Brown . ltogcat - new mapping 1757f7efa3c7SJed Brown 17589d90f715SBarry Smith Note: this currently always returns a mapping with block size of 1 17599d90f715SBarry Smith 17609d90f715SBarry Smith Developer Note: If all the input mapping have the same block size we could easily handle that as a special case 17619d90f715SBarry Smith 1762f7efa3c7SJed Brown Level: advanced 1763f7efa3c7SJed Brown 1764f7efa3c7SJed Brown .seealso: ISLocalToGlobalMappingCreate() 1765f7efa3c7SJed Brown @*/ 1766f7efa3c7SJed Brown PetscErrorCode ISLocalToGlobalMappingConcatenate(MPI_Comm comm,PetscInt n,const ISLocalToGlobalMapping ltogs[],ISLocalToGlobalMapping *ltogcat) 1767f7efa3c7SJed Brown { 1768f7efa3c7SJed Brown PetscInt i,cnt,m,*idx; 1769f7efa3c7SJed Brown PetscErrorCode ierr; 1770f7efa3c7SJed Brown 1771f7efa3c7SJed Brown PetscFunctionBegin; 1772f7efa3c7SJed Brown if (n < 0) SETERRQ1(comm,PETSC_ERR_ARG_OUTOFRANGE,"Must have a non-negative number of mappings, given %D",n); 1773f7efa3c7SJed Brown if (n > 0) PetscValidPointer(ltogs,3); 1774f7efa3c7SJed Brown for (i=0; i<n; i++) PetscValidHeaderSpecific(ltogs[i],IS_LTOGM_CLASSID,3); 1775f7efa3c7SJed Brown PetscValidPointer(ltogcat,4); 1776f7efa3c7SJed Brown for (cnt=0,i=0; i<n; i++) { 1777f7efa3c7SJed Brown ierr = ISLocalToGlobalMappingGetSize(ltogs[i],&m);CHKERRQ(ierr); 1778f7efa3c7SJed Brown cnt += m; 1779f7efa3c7SJed Brown } 1780785e854fSJed Brown ierr = PetscMalloc1(cnt,&idx);CHKERRQ(ierr); 1781f7efa3c7SJed Brown for (cnt=0,i=0; i<n; i++) { 1782f7efa3c7SJed Brown const PetscInt *subidx; 1783f7efa3c7SJed Brown ierr = ISLocalToGlobalMappingGetSize(ltogs[i],&m);CHKERRQ(ierr); 1784f7efa3c7SJed Brown ierr = ISLocalToGlobalMappingGetIndices(ltogs[i],&subidx);CHKERRQ(ierr); 1785f7efa3c7SJed Brown ierr = PetscMemcpy(&idx[cnt],subidx,m*sizeof(PetscInt));CHKERRQ(ierr); 1786f7efa3c7SJed Brown ierr = ISLocalToGlobalMappingRestoreIndices(ltogs[i],&subidx);CHKERRQ(ierr); 1787f7efa3c7SJed Brown cnt += m; 1788f7efa3c7SJed Brown } 1789f0413b6fSBarry Smith ierr = ISLocalToGlobalMappingCreate(comm,1,cnt,idx,PETSC_OWN_POINTER,ltogcat);CHKERRQ(ierr); 1790f7efa3c7SJed Brown PetscFunctionReturn(0); 1791f7efa3c7SJed Brown } 179204a59952SBarry Smith 1793413f72f0SBarry Smith /*MC 1794413f72f0SBarry Smith ISLOCALTOGLOBALMAPPINGBASIC - basic implementation of the ISLocalToGlobalMapping object. When ISGlobalToLocalMappingApply() is 1795413f72f0SBarry Smith used this is good for only small and moderate size problems. 1796413f72f0SBarry Smith 1797413f72f0SBarry Smith Options Database Keys: 1798*a2b725a8SWilliam Gropp . -islocaltoglobalmapping_type basic - select this method 1799413f72f0SBarry Smith 1800413f72f0SBarry Smith Level: beginner 1801413f72f0SBarry Smith 1802413f72f0SBarry Smith .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingSetType(), ISLOCALTOGLOBALMAPPINGHASH 1803413f72f0SBarry Smith M*/ 1804413f72f0SBarry Smith PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate_Basic(ISLocalToGlobalMapping ltog) 1805413f72f0SBarry Smith { 1806413f72f0SBarry Smith PetscFunctionBegin; 1807413f72f0SBarry Smith ltog->ops->globaltolocalmappingapply = ISGlobalToLocalMappingApply_Basic; 1808413f72f0SBarry Smith ltog->ops->globaltolocalmappingsetup = ISGlobalToLocalMappingSetUp_Basic; 1809413f72f0SBarry Smith ltog->ops->globaltolocalmappingapplyblock = ISGlobalToLocalMappingApplyBlock_Basic; 1810413f72f0SBarry Smith ltog->ops->destroy = ISLocalToGlobalMappingDestroy_Basic; 1811413f72f0SBarry Smith PetscFunctionReturn(0); 1812413f72f0SBarry Smith } 1813413f72f0SBarry Smith 1814413f72f0SBarry Smith /*MC 1815413f72f0SBarry Smith ISLOCALTOGLOBALMAPPINGHASH - hash implementation of the ISLocalToGlobalMapping object. When ISGlobalToLocalMappingApply() is 1816ed56e8eaSBarry Smith used this is good for large memory problems. 1817413f72f0SBarry Smith 1818413f72f0SBarry Smith Options Database Keys: 1819*a2b725a8SWilliam Gropp . -islocaltoglobalmapping_type hash - select this method 1820413f72f0SBarry Smith 1821ed56e8eaSBarry Smith 182295452b02SPatrick Sanan Notes: 182395452b02SPatrick Sanan This is selected automatically for large problems if the user does not set the type. 1824ed56e8eaSBarry Smith 1825413f72f0SBarry Smith Level: beginner 1826413f72f0SBarry Smith 1827413f72f0SBarry Smith .seealso: ISLocalToGlobalMappingCreate(), ISLocalToGlobalMappingSetType(), ISLOCALTOGLOBALMAPPINGHASH 1828413f72f0SBarry Smith M*/ 1829413f72f0SBarry Smith PETSC_EXTERN PetscErrorCode ISLocalToGlobalMappingCreate_Hash(ISLocalToGlobalMapping ltog) 1830413f72f0SBarry Smith { 1831413f72f0SBarry Smith PetscFunctionBegin; 1832413f72f0SBarry Smith ltog->ops->globaltolocalmappingapply = ISGlobalToLocalMappingApply_Hash; 1833413f72f0SBarry Smith ltog->ops->globaltolocalmappingsetup = ISGlobalToLocalMappingSetUp_Hash; 1834413f72f0SBarry Smith ltog->ops->globaltolocalmappingapplyblock = ISGlobalToLocalMappingApplyBlock_Hash; 1835413f72f0SBarry Smith ltog->ops->destroy = ISLocalToGlobalMappingDestroy_Hash; 1836413f72f0SBarry Smith PetscFunctionReturn(0); 1837413f72f0SBarry Smith } 1838413f72f0SBarry Smith 1839413f72f0SBarry Smith 1840413f72f0SBarry Smith /*@C 1841413f72f0SBarry Smith ISLocalToGlobalMappingRegister - Adds a method for applying a global to local mapping with an ISLocalToGlobalMapping 1842413f72f0SBarry Smith 1843413f72f0SBarry Smith Not Collective 1844413f72f0SBarry Smith 1845413f72f0SBarry Smith Input Parameters: 1846413f72f0SBarry Smith + sname - name of a new method 1847413f72f0SBarry Smith - routine_create - routine to create method context 1848413f72f0SBarry Smith 1849413f72f0SBarry Smith Notes: 1850ed56e8eaSBarry Smith ISLocalToGlobalMappingRegister() may be called multiple times to add several user-defined mappings. 1851413f72f0SBarry Smith 1852413f72f0SBarry Smith Sample usage: 1853413f72f0SBarry Smith .vb 1854413f72f0SBarry Smith ISLocalToGlobalMappingRegister("my_mapper",MyCreate); 1855413f72f0SBarry Smith .ve 1856413f72f0SBarry Smith 1857ed56e8eaSBarry Smith Then, your mapping can be chosen with the procedural interface via 1858413f72f0SBarry Smith $ ISLocalToGlobalMappingSetType(ltog,"my_mapper") 1859413f72f0SBarry Smith or at runtime via the option 1860ed56e8eaSBarry Smith $ -islocaltoglobalmapping_type my_mapper 1861413f72f0SBarry Smith 1862413f72f0SBarry Smith Level: advanced 1863413f72f0SBarry Smith 1864413f72f0SBarry Smith .keywords: ISLocalToGlobalMappingType, register 1865413f72f0SBarry Smith 1866413f72f0SBarry Smith .seealso: ISLocalToGlobalMappingRegisterAll(), ISLocalToGlobalMappingRegisterDestroy(), ISLOCALTOGLOBALMAPPINGBASIC, ISLOCALTOGLOBALMAPPINGHASH 1867413f72f0SBarry Smith 1868413f72f0SBarry Smith @*/ 1869413f72f0SBarry Smith PetscErrorCode ISLocalToGlobalMappingRegister(const char sname[],PetscErrorCode (*function)(ISLocalToGlobalMapping)) 1870413f72f0SBarry Smith { 1871413f72f0SBarry Smith PetscErrorCode ierr; 1872413f72f0SBarry Smith 1873413f72f0SBarry Smith PetscFunctionBegin; 18741d36bdfdSBarry Smith ierr = ISInitializePackage();CHKERRQ(ierr); 1875413f72f0SBarry Smith ierr = PetscFunctionListAdd(&ISLocalToGlobalMappingList,sname,function);CHKERRQ(ierr); 1876413f72f0SBarry Smith PetscFunctionReturn(0); 1877413f72f0SBarry Smith } 1878413f72f0SBarry Smith 1879413f72f0SBarry Smith /*@C 1880ed56e8eaSBarry Smith ISLocalToGlobalMappingSetType - Builds ISLocalToGlobalMapping for a particular global to local mapping approach. 1881413f72f0SBarry Smith 1882413f72f0SBarry Smith Logically Collective on ISLocalToGlobalMapping 1883413f72f0SBarry Smith 1884413f72f0SBarry Smith Input Parameters: 1885413f72f0SBarry Smith + ltog - the ISLocalToGlobalMapping object 1886413f72f0SBarry Smith - type - a known method 1887413f72f0SBarry Smith 1888413f72f0SBarry Smith Options Database Key: 1889ed56e8eaSBarry Smith . -islocaltoglobalmapping_type <method> - Sets the method; use -help for a list 1890413f72f0SBarry Smith of available methods (for instance, basic or hash) 1891413f72f0SBarry Smith 1892413f72f0SBarry Smith Notes: 1893413f72f0SBarry Smith See "petsc/include/petscis.h" for available methods 1894413f72f0SBarry Smith 1895413f72f0SBarry Smith Normally, it is best to use the ISLocalToGlobalMappingSetFromOptions() command and 1896413f72f0SBarry Smith then set the ISLocalToGlobalMapping type from the options database rather than by using 1897413f72f0SBarry Smith this routine. 1898413f72f0SBarry Smith 1899413f72f0SBarry Smith Level: intermediate 1900413f72f0SBarry Smith 1901413f72f0SBarry Smith Developer Note: ISLocalToGlobalMappingRegister() is used to add new types to ISLocalToGlobalMappingList from which they 1902413f72f0SBarry Smith are accessed by ISLocalToGlobalMappingSetType(). 1903413f72f0SBarry Smith 1904413f72f0SBarry Smith .keywords: ISLocalToGlobalMapping, set, method 1905413f72f0SBarry Smith 1906413f72f0SBarry Smith .seealso: PCSetType(), ISLocalToGlobalMappingType, ISLocalToGlobalMappingRegister(), ISLocalToGlobalMappingCreate() 1907413f72f0SBarry Smith 1908413f72f0SBarry Smith @*/ 1909413f72f0SBarry Smith PetscErrorCode ISLocalToGlobalMappingSetType(ISLocalToGlobalMapping ltog, ISLocalToGlobalMappingType type) 1910413f72f0SBarry Smith { 1911413f72f0SBarry Smith PetscErrorCode ierr,(*r)(ISLocalToGlobalMapping); 1912413f72f0SBarry Smith PetscBool match; 1913413f72f0SBarry Smith 1914413f72f0SBarry Smith PetscFunctionBegin; 1915413f72f0SBarry Smith PetscValidHeaderSpecific(ltog,IS_LTOGM_CLASSID,1); 1916413f72f0SBarry Smith PetscValidCharPointer(type,2); 1917413f72f0SBarry Smith 1918413f72f0SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)ltog,type,&match);CHKERRQ(ierr); 1919413f72f0SBarry Smith if (match) PetscFunctionReturn(0); 1920413f72f0SBarry Smith 1921413f72f0SBarry Smith ierr = PetscFunctionListFind(ISLocalToGlobalMappingList,type,&r);CHKERRQ(ierr); 1922413f72f0SBarry Smith if (!r) SETERRQ1(PetscObjectComm((PetscObject)ltog),PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested ISLocalToGlobalMapping type %s",type); 1923413f72f0SBarry Smith /* Destroy the previous private LTOG context */ 1924413f72f0SBarry Smith if (ltog->ops->destroy) { 1925413f72f0SBarry Smith ierr = (*ltog->ops->destroy)(ltog);CHKERRQ(ierr); 1926413f72f0SBarry Smith ltog->ops->destroy = NULL; 1927413f72f0SBarry Smith } 1928413f72f0SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)ltog,type);CHKERRQ(ierr); 1929413f72f0SBarry Smith ierr = (*r)(ltog);CHKERRQ(ierr); 1930413f72f0SBarry Smith PetscFunctionReturn(0); 1931413f72f0SBarry Smith } 1932413f72f0SBarry Smith 1933413f72f0SBarry Smith PetscBool ISLocalToGlobalMappingRegisterAllCalled = PETSC_FALSE; 1934413f72f0SBarry Smith 1935413f72f0SBarry Smith /*@C 1936413f72f0SBarry Smith ISLocalToGlobalMappingRegisterAll - Registers all of the local to global mapping components in the IS package. 1937413f72f0SBarry Smith 1938413f72f0SBarry Smith Not Collective 1939413f72f0SBarry Smith 1940413f72f0SBarry Smith Level: advanced 1941413f72f0SBarry Smith 1942413f72f0SBarry Smith .keywords: IS, register, all 1943413f72f0SBarry Smith .seealso: ISRegister(), ISLocalToGlobalRegister() 1944413f72f0SBarry Smith @*/ 1945413f72f0SBarry Smith PetscErrorCode ISLocalToGlobalMappingRegisterAll(void) 1946413f72f0SBarry Smith { 1947413f72f0SBarry Smith PetscErrorCode ierr; 1948413f72f0SBarry Smith 1949413f72f0SBarry Smith PetscFunctionBegin; 1950413f72f0SBarry Smith if (ISLocalToGlobalMappingRegisterAllCalled) PetscFunctionReturn(0); 1951413f72f0SBarry Smith ISLocalToGlobalMappingRegisterAllCalled = PETSC_TRUE; 1952413f72f0SBarry Smith 1953413f72f0SBarry Smith ierr = ISLocalToGlobalMappingRegister(ISLOCALTOGLOBALMAPPINGBASIC, ISLocalToGlobalMappingCreate_Basic);CHKERRQ(ierr); 1954413f72f0SBarry Smith ierr = ISLocalToGlobalMappingRegister(ISLOCALTOGLOBALMAPPINGHASH, ISLocalToGlobalMappingCreate_Hash);CHKERRQ(ierr); 1955413f72f0SBarry Smith PetscFunctionReturn(0); 1956413f72f0SBarry Smith } 195704a59952SBarry Smith 1958