1c6db04a5SJed Brown #include <petscdmsliced.h> /*I "petscdmsliced.h" I*/ 207475bc1SBarry Smith #include <petscmat.h> 3af0996ceSBarry Smith #include <petsc/private/dmimpl.h> 447c6ae99SBarry Smith 547c6ae99SBarry Smith /* CSR storage of the nonzero structure of a bs*bs matrix */ 647c6ae99SBarry Smith typedef struct { 747c6ae99SBarry Smith PetscInt bs, nz, *i, *j; 80c010503SBarry Smith } DMSlicedBlockFills; 947c6ae99SBarry Smith 1047c6ae99SBarry Smith typedef struct { 1147c6ae99SBarry Smith PetscInt bs, n, N, Nghosts, *ghosts; 1247c6ae99SBarry Smith PetscInt d_nz, o_nz, *d_nnz, *o_nnz; 130c010503SBarry Smith DMSlicedBlockFills *dfill, *ofill; 1447c6ae99SBarry Smith } DM_Sliced; 1547c6ae99SBarry Smith 16*a4e35b19SJacob Faibussowitsch static PetscErrorCode DMCreateMatrix_Sliced(DM dm, Mat *J) 17d71ae5a4SJacob Faibussowitsch { 1847c6ae99SBarry Smith PetscInt *globals, *sd_nnz, *so_nnz, rstart, bs, i; 1945b6f7e9SBarry Smith ISLocalToGlobalMapping lmap; 20f77a5242SKarl Rupp void (*aij)(void) = NULL; 2147c6ae99SBarry Smith DM_Sliced *slice = (DM_Sliced *)dm->data; 2247c6ae99SBarry Smith 2347c6ae99SBarry Smith PetscFunctionBegin; 2447c6ae99SBarry Smith bs = slice->bs; 259566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)dm), J)); 269566063dSJacob Faibussowitsch PetscCall(MatSetSizes(*J, slice->n * bs, slice->n * bs, PETSC_DETERMINE, PETSC_DETERMINE)); 279566063dSJacob Faibussowitsch PetscCall(MatSetBlockSize(*J, bs)); 289566063dSJacob Faibussowitsch PetscCall(MatSetType(*J, dm->mattype)); 299566063dSJacob Faibussowitsch PetscCall(MatSeqBAIJSetPreallocation(*J, bs, slice->d_nz, slice->d_nnz)); 309566063dSJacob Faibussowitsch PetscCall(MatMPIBAIJSetPreallocation(*J, bs, slice->d_nz, slice->d_nnz, slice->o_nz, slice->o_nnz)); 3147c6ae99SBarry Smith /* In general, we have to do extra work to preallocate for scalar (AIJ) matrices so we check whether it will do any 3247c6ae99SBarry Smith * good before going on with it. */ 339566063dSJacob Faibussowitsch PetscCall(PetscObjectQueryFunction((PetscObject)*J, "MatMPIAIJSetPreallocation_C", &aij)); 3448a46eb9SPierre Jolivet if (!aij) PetscCall(PetscObjectQueryFunction((PetscObject)*J, "MatSeqAIJSetPreallocation_C", &aij)); 3547c6ae99SBarry Smith if (aij) { 3647c6ae99SBarry Smith if (bs == 1) { 379566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(*J, slice->d_nz, slice->d_nnz)); 389566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(*J, slice->d_nz, slice->d_nnz, slice->o_nz, slice->o_nnz)); 3947c6ae99SBarry Smith } else if (!slice->d_nnz) { 409566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(*J, slice->d_nz * bs, NULL)); 419566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(*J, slice->d_nz * bs, NULL, slice->o_nz * bs, NULL)); 4247c6ae99SBarry Smith } else { 430c010503SBarry Smith /* The user has provided preallocation per block-row, convert it to per scalar-row respecting DMSlicedSetBlockFills() if applicable */ 449566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(slice->n * bs, &sd_nnz, (!!slice->o_nnz) * slice->n * bs, &so_nnz)); 4547c6ae99SBarry Smith for (i = 0; i < slice->n * bs; i++) { 469371c9d4SSatish Balay sd_nnz[i] = (slice->d_nnz[i / bs] - 1) * (slice->ofill ? slice->ofill->i[i % bs + 1] - slice->ofill->i[i % bs] : bs) + (slice->dfill ? slice->dfill->i[i % bs + 1] - slice->dfill->i[i % bs] : bs); 47ad540459SPierre Jolivet if (so_nnz) so_nnz[i] = slice->o_nnz[i / bs] * (slice->ofill ? slice->ofill->i[i % bs + 1] - slice->ofill->i[i % bs] : bs); 4847c6ae99SBarry Smith } 499566063dSJacob Faibussowitsch PetscCall(MatSeqAIJSetPreallocation(*J, slice->d_nz * bs, sd_nnz)); 509566063dSJacob Faibussowitsch PetscCall(MatMPIAIJSetPreallocation(*J, slice->d_nz * bs, sd_nnz, slice->o_nz * bs, so_nnz)); 519566063dSJacob Faibussowitsch PetscCall(PetscFree2(sd_nnz, so_nnz)); 5247c6ae99SBarry Smith } 5347c6ae99SBarry Smith } 5447c6ae99SBarry Smith 5547c6ae99SBarry Smith /* Set up the local to global map. For the scalar map, we have to translate to entry-wise indexing instead of block-wise. */ 569566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(slice->n + slice->Nghosts, &globals)); 579566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(*J, &rstart, NULL)); 5845b6f7e9SBarry Smith for (i = 0; i < slice->n; i++) globals[i] = rstart / bs + i; 598865f1eaSKarl Rupp 60ad540459SPierre Jolivet for (i = 0; i < slice->Nghosts; i++) globals[slice->n + i] = slice->ghosts[i]; 619566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreate(PETSC_COMM_SELF, bs, slice->n + slice->Nghosts, globals, PETSC_OWN_POINTER, &lmap)); 629566063dSJacob Faibussowitsch PetscCall(MatSetLocalToGlobalMapping(*J, lmap, lmap)); 639566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(&lmap)); 649566063dSJacob Faibussowitsch PetscCall(MatSetDM(*J, dm)); 653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6647c6ae99SBarry Smith } 6747c6ae99SBarry Smith 6847c6ae99SBarry Smith /*@C 690c010503SBarry Smith DMSlicedSetGhosts - Sets the global indices of other processes elements that will 7047c6ae99SBarry Smith be ghosts on this process 7147c6ae99SBarry Smith 7247c6ae99SBarry Smith Not Collective 7347c6ae99SBarry Smith 7447c6ae99SBarry Smith Input Parameters: 7560225df5SJacob Faibussowitsch + dm - the `DMSLICED` object 7647c6ae99SBarry Smith . bs - block size 7747c6ae99SBarry Smith . nlocal - number of local (owned, non-ghost) blocks 7847c6ae99SBarry Smith . Nghosts - number of ghost blocks on this process 7947c6ae99SBarry Smith - ghosts - global indices of each ghost block 8047c6ae99SBarry Smith 8147c6ae99SBarry Smith Level: advanced 8247c6ae99SBarry Smith 832fe279fdSBarry Smith .seealso: `DM`, `DMSLICED`, `DMDestroy()`, `DMCreateGlobalVector()` 8447c6ae99SBarry Smith @*/ 85d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSlicedSetGhosts(DM dm, PetscInt bs, PetscInt nlocal, PetscInt Nghosts, const PetscInt ghosts[]) 86d71ae5a4SJacob Faibussowitsch { 8747c6ae99SBarry Smith DM_Sliced *slice = (DM_Sliced *)dm->data; 8847c6ae99SBarry Smith 8947c6ae99SBarry Smith PetscFunctionBegin; 900c010503SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 919566063dSJacob Faibussowitsch PetscCall(PetscFree(slice->ghosts)); 929566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Nghosts, &slice->ghosts)); 939566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(slice->ghosts, ghosts, Nghosts)); 9447c6ae99SBarry Smith slice->bs = bs; 9547c6ae99SBarry Smith slice->n = nlocal; 9647c6ae99SBarry Smith slice->Nghosts = Nghosts; 973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9847c6ae99SBarry Smith } 9947c6ae99SBarry Smith 10047c6ae99SBarry Smith /*@C 1012fe279fdSBarry Smith DMSlicedSetPreallocation - sets the matrix memory preallocation for matrices computed by `DMSLICED` 10247c6ae99SBarry Smith 10347c6ae99SBarry Smith Not Collective 10447c6ae99SBarry Smith 10547c6ae99SBarry Smith Input Parameters: 10660225df5SJacob Faibussowitsch + dm - the `DM` object 10747c6ae99SBarry Smith . d_nz - number of block nonzeros per block row in diagonal portion of local 10847c6ae99SBarry Smith submatrix (same for all local rows) 10947c6ae99SBarry Smith . d_nnz - array containing the number of block nonzeros in the various block rows 11047c6ae99SBarry Smith of the in diagonal portion of the local (possibly different for each block 1112fe279fdSBarry Smith row) or `NULL`. 11247c6ae99SBarry Smith . o_nz - number of block nonzeros per block row in the off-diagonal portion of local 11347c6ae99SBarry Smith submatrix (same for all local rows). 11447c6ae99SBarry Smith - o_nnz - array containing the number of nonzeros in the various block rows of the 11547c6ae99SBarry Smith off-diagonal portion of the local submatrix (possibly different for 11620f4b53cSBarry Smith each block row) or `NULL`. 11747c6ae99SBarry Smith 11847c6ae99SBarry Smith Level: advanced 11947c6ae99SBarry Smith 12020f4b53cSBarry Smith Note: 1212fe279fdSBarry Smith See `MatMPIBAIJSetPreallocation()` for more details on preallocation. If a scalar matrix (`MATAIJ`) is 12220f4b53cSBarry Smith obtained with `DMSlicedGetMatrix()`, the correct preallocation will be set, respecting `DMSlicedSetBlockFills()`. 12320f4b53cSBarry Smith 1242fe279fdSBarry Smith .seealso: `DM`, `DMSLICED`, `DMDestroy()`, `DMCreateGlobalVector()`, `MatMPIAIJSetPreallocation()`, 125db781477SPatrick Sanan `MatMPIBAIJSetPreallocation()`, `DMSlicedGetMatrix()`, `DMSlicedSetBlockFills()` 12647c6ae99SBarry Smith @*/ 127d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSlicedSetPreallocation(DM dm, PetscInt d_nz, const PetscInt d_nnz[], PetscInt o_nz, const PetscInt o_nnz[]) 128d71ae5a4SJacob Faibussowitsch { 12947c6ae99SBarry Smith DM_Sliced *slice = (DM_Sliced *)dm->data; 13047c6ae99SBarry Smith 13147c6ae99SBarry Smith PetscFunctionBegin; 13247c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 13347c6ae99SBarry Smith slice->d_nz = d_nz; 13447c6ae99SBarry Smith slice->d_nnz = (PetscInt *)d_nnz; 13547c6ae99SBarry Smith slice->o_nz = o_nz; 13647c6ae99SBarry Smith slice->o_nnz = (PetscInt *)o_nnz; 1373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13847c6ae99SBarry Smith } 13947c6ae99SBarry Smith 140d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSlicedSetBlockFills_Private(PetscInt bs, const PetscInt *fill, DMSlicedBlockFills **inf) 141d71ae5a4SJacob Faibussowitsch { 14247c6ae99SBarry Smith PetscInt i, j, nz, *fi, *fj; 1430c010503SBarry Smith DMSlicedBlockFills *f; 14447c6ae99SBarry Smith 14547c6ae99SBarry Smith PetscFunctionBegin; 1464f572ea9SToby Isaac PetscAssertPointer(inf, 3); 1479566063dSJacob Faibussowitsch if (*inf) PetscCall(PetscFree3(*inf, (*inf)->i, (*inf)->j)); 1483ba16761SJacob Faibussowitsch if (!fill) PetscFunctionReturn(PETSC_SUCCESS); 1499371c9d4SSatish Balay for (i = 0, nz = 0; i < bs * bs; i++) 1509371c9d4SSatish Balay if (fill[i]) nz++; 1519566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(1, &f, bs + 1, &fi, nz, &fj)); 15247c6ae99SBarry Smith f->bs = bs; 15347c6ae99SBarry Smith f->nz = nz; 15447c6ae99SBarry Smith f->i = fi; 15547c6ae99SBarry Smith f->j = fj; 15647c6ae99SBarry Smith for (i = 0, nz = 0; i < bs; i++) { 15747c6ae99SBarry Smith fi[i] = nz; 1589371c9d4SSatish Balay for (j = 0; j < bs; j++) 1599371c9d4SSatish Balay if (fill[i * bs + j]) fj[nz++] = j; 16047c6ae99SBarry Smith } 16147c6ae99SBarry Smith fi[i] = nz; 16247c6ae99SBarry Smith *inf = f; 1633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16447c6ae99SBarry Smith } 16547c6ae99SBarry Smith 16647c6ae99SBarry Smith /*@C 1670c010503SBarry Smith DMSlicedSetBlockFills - Sets the fill pattern in each block for a multi-component problem 16820f4b53cSBarry Smith of the matrix returned by `DMSlicedGetMatrix()`. 16947c6ae99SBarry Smith 17020f4b53cSBarry Smith Logically Collective 17147c6ae99SBarry Smith 172d8d19677SJose E. Roman Input Parameters: 17360225df5SJacob Faibussowitsch + dm - the `DM` object 17420f4b53cSBarry Smith . dfill - the fill pattern in the diagonal block (may be `NULL`, means use dense block) 17547c6ae99SBarry Smith - ofill - the fill pattern in the off-diagonal blocks 17647c6ae99SBarry Smith 17747c6ae99SBarry Smith Level: advanced 17847c6ae99SBarry Smith 17920f4b53cSBarry Smith Note: 18020f4b53cSBarry Smith This only makes sense for multicomponent problems using scalar matrix formats (AIJ). 18120f4b53cSBarry Smith See `DMDASetBlockFills()` for example usage. 18220f4b53cSBarry Smith 1832fe279fdSBarry Smith .seealso: `DM`, `DMSLICED`, `DMSlicedGetMatrix()`, `DMDASetBlockFills()` 18447c6ae99SBarry Smith @*/ 185d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSlicedSetBlockFills(DM dm, const PetscInt *dfill, const PetscInt *ofill) 186d71ae5a4SJacob Faibussowitsch { 18747c6ae99SBarry Smith DM_Sliced *slice = (DM_Sliced *)dm->data; 18847c6ae99SBarry Smith 18947c6ae99SBarry Smith PetscFunctionBegin; 1900c010503SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1919566063dSJacob Faibussowitsch PetscCall(DMSlicedSetBlockFills_Private(slice->bs, dfill, &slice->dfill)); 1929566063dSJacob Faibussowitsch PetscCall(DMSlicedSetBlockFills_Private(slice->bs, ofill, &slice->ofill)); 1933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 19447c6ae99SBarry Smith } 19547c6ae99SBarry Smith 196d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDestroy_Sliced(DM dm) 197d71ae5a4SJacob Faibussowitsch { 19847c6ae99SBarry Smith DM_Sliced *slice = (DM_Sliced *)dm->data; 19947c6ae99SBarry Smith 20047c6ae99SBarry Smith PetscFunctionBegin; 2019566063dSJacob Faibussowitsch PetscCall(PetscFree(slice->ghosts)); 2029566063dSJacob Faibussowitsch if (slice->dfill) PetscCall(PetscFree3(slice->dfill, slice->dfill->i, slice->dfill->j)); 2039566063dSJacob Faibussowitsch if (slice->ofill) PetscCall(PetscFree3(slice->ofill, slice->ofill->i, slice->ofill->j)); 204435a35e8SMatthew G Knepley /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */ 2059566063dSJacob Faibussowitsch PetscCall(PetscFree(slice)); 2063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 20747c6ae99SBarry Smith } 20847c6ae99SBarry Smith 209d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCreateGlobalVector_Sliced(DM dm, Vec *gvec) 210d71ae5a4SJacob Faibussowitsch { 21147c6ae99SBarry Smith DM_Sliced *slice = (DM_Sliced *)dm->data; 21247c6ae99SBarry Smith 21347c6ae99SBarry Smith PetscFunctionBegin; 21447c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2154f572ea9SToby Isaac PetscAssertPointer(gvec, 2); 216ea78f98cSLisandro Dalcin *gvec = NULL; 2179566063dSJacob Faibussowitsch PetscCall(VecCreateGhostBlock(PetscObjectComm((PetscObject)dm), slice->bs, slice->n * slice->bs, PETSC_DETERMINE, slice->Nghosts, slice->ghosts, gvec)); 2189566063dSJacob Faibussowitsch PetscCall(VecSetDM(*gvec, dm)); 2193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 22047c6ae99SBarry Smith } 22147c6ae99SBarry Smith 222d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMGlobalToLocalBegin_Sliced(DM da, Vec g, InsertMode mode, Vec l) 223d71ae5a4SJacob Faibussowitsch { 2243efe6655SBarry Smith PetscBool flg; 2253efe6655SBarry Smith 2263efe6655SBarry Smith PetscFunctionBegin; 2279566063dSJacob Faibussowitsch PetscCall(VecGhostIsLocalForm(g, l, &flg)); 22828b400f6SJacob Faibussowitsch PetscCheck(flg, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONG, "Local vector is not local form of global vector"); 2299566063dSJacob Faibussowitsch PetscCall(VecGhostUpdateEnd(g, mode, SCATTER_FORWARD)); 2309566063dSJacob Faibussowitsch PetscCall(VecGhostUpdateBegin(g, mode, SCATTER_FORWARD)); 2313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2323efe6655SBarry Smith } 2333efe6655SBarry Smith 234d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMGlobalToLocalEnd_Sliced(DM da, Vec g, InsertMode mode, Vec l) 235d71ae5a4SJacob Faibussowitsch { 2363efe6655SBarry Smith PetscBool flg; 2373efe6655SBarry Smith 2383efe6655SBarry Smith PetscFunctionBegin; 2399566063dSJacob Faibussowitsch PetscCall(VecGhostIsLocalForm(g, l, &flg)); 24028b400f6SJacob Faibussowitsch PetscCheck(flg, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONG, "Local vector is not local form of global vector"); 2419566063dSJacob Faibussowitsch PetscCall(VecGhostUpdateEnd(g, mode, SCATTER_FORWARD)); 2423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2433efe6655SBarry Smith } 2443efe6655SBarry Smith 2453e0de288SBarry Smith /*MC 2462fe279fdSBarry Smith DMSLICED = "sliced" - A `DM` object that is used to manage data for a general graph. Uses `VecCreateGhost()` ghosted vectors for storing the fields 2473e0de288SBarry Smith 2482fe279fdSBarry Smith See `DMCreateSliced()` for details. 2493e0de288SBarry Smith 2503e0de288SBarry Smith Level: intermediate 2513e0de288SBarry Smith 252db781477SPatrick Sanan .seealso: `DMType`, `DMCOMPOSITE`, `DMCreateSliced()`, `DMCreate()` 2533e0de288SBarry Smith M*/ 2543e0de288SBarry Smith 255d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMCreate_Sliced(DM p) 256d71ae5a4SJacob Faibussowitsch { 257a4121054SBarry Smith DM_Sliced *slice; 258a4121054SBarry Smith 259a4121054SBarry Smith PetscFunctionBegin; 2604dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&slice)); 261a4121054SBarry Smith p->data = slice; 262a4121054SBarry Smith 263a4121054SBarry Smith p->ops->createglobalvector = DMCreateGlobalVector_Sliced; 26425296bd5SBarry Smith p->ops->creatematrix = DMCreateMatrix_Sliced; 2653efe6655SBarry Smith p->ops->globaltolocalbegin = DMGlobalToLocalBegin_Sliced; 2663efe6655SBarry Smith p->ops->globaltolocalend = DMGlobalToLocalEnd_Sliced; 267a4121054SBarry Smith p->ops->destroy = DMDestroy_Sliced; 2683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 269a4121054SBarry Smith } 270a4121054SBarry Smith 27147c6ae99SBarry Smith /*@C 27220f4b53cSBarry Smith DMSlicedCreate - Creates a `DM` object, used to manage data for a unstructured problem 27347c6ae99SBarry Smith 274d083f849SBarry Smith Collective 2750c010503SBarry Smith 276d8d19677SJose E. Roman Input Parameters: 2773efe6655SBarry Smith + comm - the processors that will share the global vector 2783efe6655SBarry Smith . bs - the block size 2793efe6655SBarry Smith . nlocal - number of vector entries on this process 2803efe6655SBarry Smith . Nghosts - number of ghost points needed on this process 2813efe6655SBarry Smith . ghosts - global indices of all ghost points for this process 2823efe6655SBarry Smith . d_nnz - matrix preallocation information representing coupling within this process 2833efe6655SBarry Smith - o_nnz - matrix preallocation information representing coupling between this process and other processes 2840c010503SBarry Smith 2852fe279fdSBarry Smith Output Parameter: 28660225df5SJacob Faibussowitsch . dm - the slice object 2870c010503SBarry Smith 28820f4b53cSBarry Smith Level: advanced 28920f4b53cSBarry Smith 2903efe6655SBarry Smith Notes: 29120f4b53cSBarry Smith This `DM` does not support `DMCreateLocalVector()`, `DMGlobalToLocalBegin()`, and `DMGlobalToLocalEnd()` instead one directly uses 29220f4b53cSBarry Smith `VecGhostGetLocalForm()` and `VecGhostRestoreLocalForm()` to access the local representation and `VecGhostUpdateBegin()` and `VecGhostUpdateEnd()` to update 2933efe6655SBarry Smith the ghost points. 2943efe6655SBarry Smith 29520f4b53cSBarry Smith One can use `DMGlobalToLocalBegin()`, and `DMGlobalToLocalEnd()` instead of `VecGhostUpdateBegin()` and `VecGhostUpdateEnd()`. 2960c010503SBarry Smith 29760225df5SJacob Faibussowitsch .seealso: `DM`, `DMSLICED`, `DMDestroy()`, `DMCreateGlobalVector()`, `DMSetType()`, `DMSlicedSetGhosts()`, `DMSlicedSetPreallocation()`, 2982fe279fdSBarry Smith `VecGhostUpdateBegin()`, `VecGhostUpdateEnd()`, 299db781477SPatrick Sanan `VecGhostGetLocalForm()`, `VecGhostRestoreLocalForm()` 3000c010503SBarry Smith @*/ 301d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSlicedCreate(MPI_Comm comm, PetscInt bs, PetscInt nlocal, PetscInt Nghosts, const PetscInt ghosts[], const PetscInt d_nnz[], const PetscInt o_nnz[], DM *dm) 302d71ae5a4SJacob Faibussowitsch { 3030c010503SBarry Smith PetscFunctionBegin; 3044f572ea9SToby Isaac PetscAssertPointer(dm, 8); 3059566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 3069566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMSLICED)); 3079566063dSJacob Faibussowitsch PetscCall(DMSlicedSetGhosts(*dm, bs, nlocal, Nghosts, ghosts)); 3081baa6e33SBarry Smith if (d_nnz) PetscCall(DMSlicedSetPreallocation(*dm, 0, d_nnz, 0, o_nnz)); 3093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3100c010503SBarry Smith } 311