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*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMCreateMatrix_Sliced(DM dm, Mat *J) 17*d71ae5a4SJacob 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)); 6547c6ae99SBarry Smith PetscFunctionReturn(0); 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: 750c010503SBarry Smith + slice - the DM 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 83db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCreateGlobalVector()` 8447c6ae99SBarry Smith 8547c6ae99SBarry Smith @*/ 86*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSlicedSetGhosts(DM dm, PetscInt bs, PetscInt nlocal, PetscInt Nghosts, const PetscInt ghosts[]) 87*d71ae5a4SJacob Faibussowitsch { 8847c6ae99SBarry Smith DM_Sliced *slice = (DM_Sliced *)dm->data; 8947c6ae99SBarry Smith 9047c6ae99SBarry Smith PetscFunctionBegin; 910c010503SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 929566063dSJacob Faibussowitsch PetscCall(PetscFree(slice->ghosts)); 939566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Nghosts, &slice->ghosts)); 949566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(slice->ghosts, ghosts, Nghosts)); 9547c6ae99SBarry Smith slice->bs = bs; 9647c6ae99SBarry Smith slice->n = nlocal; 9747c6ae99SBarry Smith slice->Nghosts = Nghosts; 9847c6ae99SBarry Smith PetscFunctionReturn(0); 9947c6ae99SBarry Smith } 10047c6ae99SBarry Smith 10147c6ae99SBarry Smith /*@C 1020c010503SBarry Smith DMSlicedSetPreallocation - sets the matrix memory preallocation for matrices computed by DMSliced 10347c6ae99SBarry Smith 10447c6ae99SBarry Smith Not Collective 10547c6ae99SBarry Smith 10647c6ae99SBarry Smith Input Parameters: 1070c010503SBarry Smith + slice - the DM object 10847c6ae99SBarry Smith . d_nz - number of block nonzeros per block row in diagonal portion of local 10947c6ae99SBarry Smith submatrix (same for all local rows) 11047c6ae99SBarry Smith . d_nnz - array containing the number of block nonzeros in the various block rows 11147c6ae99SBarry Smith of the in diagonal portion of the local (possibly different for each block 112f77a5242SKarl Rupp row) or NULL. 11347c6ae99SBarry Smith . o_nz - number of block nonzeros per block row in the off-diagonal portion of local 11447c6ae99SBarry Smith submatrix (same for all local rows). 11547c6ae99SBarry Smith - o_nnz - array containing the number of nonzeros in the various block rows of the 11647c6ae99SBarry Smith off-diagonal portion of the local submatrix (possibly different for 117f77a5242SKarl Rupp each block row) or NULL. 11847c6ae99SBarry Smith 11947c6ae99SBarry Smith Notes: 12047c6ae99SBarry Smith See MatMPIBAIJSetPreallocation() for more details on preallocation. If a scalar matrix (AIJ) is 1210c010503SBarry Smith obtained with DMSlicedGetMatrix(), the correct preallocation will be set, respecting DMSlicedSetBlockFills(). 12247c6ae99SBarry Smith 12347c6ae99SBarry Smith Level: advanced 12447c6ae99SBarry Smith 125db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCreateGlobalVector()`, `MatMPIAIJSetPreallocation()`, 126db781477SPatrick Sanan `MatMPIBAIJSetPreallocation()`, `DMSlicedGetMatrix()`, `DMSlicedSetBlockFills()` 12747c6ae99SBarry Smith 12847c6ae99SBarry Smith @*/ 129*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSlicedSetPreallocation(DM dm, PetscInt d_nz, const PetscInt d_nnz[], PetscInt o_nz, const PetscInt o_nnz[]) 130*d71ae5a4SJacob Faibussowitsch { 13147c6ae99SBarry Smith DM_Sliced *slice = (DM_Sliced *)dm->data; 13247c6ae99SBarry Smith 13347c6ae99SBarry Smith PetscFunctionBegin; 13447c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 13547c6ae99SBarry Smith slice->d_nz = d_nz; 13647c6ae99SBarry Smith slice->d_nnz = (PetscInt *)d_nnz; 13747c6ae99SBarry Smith slice->o_nz = o_nz; 13847c6ae99SBarry Smith slice->o_nnz = (PetscInt *)o_nnz; 13947c6ae99SBarry Smith PetscFunctionReturn(0); 14047c6ae99SBarry Smith } 14147c6ae99SBarry Smith 142*d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMSlicedSetBlockFills_Private(PetscInt bs, const PetscInt *fill, DMSlicedBlockFills **inf) 143*d71ae5a4SJacob Faibussowitsch { 14447c6ae99SBarry Smith PetscInt i, j, nz, *fi, *fj; 1450c010503SBarry Smith DMSlicedBlockFills *f; 14647c6ae99SBarry Smith 14747c6ae99SBarry Smith PetscFunctionBegin; 14847c6ae99SBarry Smith PetscValidPointer(inf, 3); 1499566063dSJacob Faibussowitsch if (*inf) PetscCall(PetscFree3(*inf, (*inf)->i, (*inf)->j)); 15047c6ae99SBarry Smith if (!fill) PetscFunctionReturn(0); 1519371c9d4SSatish Balay for (i = 0, nz = 0; i < bs * bs; i++) 1529371c9d4SSatish Balay if (fill[i]) nz++; 1539566063dSJacob Faibussowitsch PetscCall(PetscMalloc3(1, &f, bs + 1, &fi, nz, &fj)); 15447c6ae99SBarry Smith f->bs = bs; 15547c6ae99SBarry Smith f->nz = nz; 15647c6ae99SBarry Smith f->i = fi; 15747c6ae99SBarry Smith f->j = fj; 15847c6ae99SBarry Smith for (i = 0, nz = 0; i < bs; i++) { 15947c6ae99SBarry Smith fi[i] = nz; 1609371c9d4SSatish Balay for (j = 0; j < bs; j++) 1619371c9d4SSatish Balay if (fill[i * bs + j]) fj[nz++] = j; 16247c6ae99SBarry Smith } 16347c6ae99SBarry Smith fi[i] = nz; 16447c6ae99SBarry Smith *inf = f; 16547c6ae99SBarry Smith PetscFunctionReturn(0); 16647c6ae99SBarry Smith } 16747c6ae99SBarry Smith 16847c6ae99SBarry Smith /*@C 1690c010503SBarry Smith DMSlicedSetBlockFills - Sets the fill pattern in each block for a multi-component problem 1700c010503SBarry Smith of the matrix returned by DMSlicedGetMatrix(). 17147c6ae99SBarry Smith 172d083f849SBarry Smith Logically Collective on dm 17347c6ae99SBarry Smith 174d8d19677SJose E. Roman Input Parameters: 1750c010503SBarry Smith + sliced - the DM object 176f77a5242SKarl Rupp . dfill - the fill pattern in the diagonal block (may be NULL, means use dense block) 17747c6ae99SBarry Smith - ofill - the fill pattern in the off-diagonal blocks 17847c6ae99SBarry Smith 17947c6ae99SBarry Smith Notes: 18047c6ae99SBarry Smith This only makes sense for multicomponent problems using scalar matrix formats (AIJ). 181aa219208SBarry Smith See DMDASetBlockFills() for example usage. 18247c6ae99SBarry Smith 18347c6ae99SBarry Smith Level: advanced 18447c6ae99SBarry Smith 185db781477SPatrick Sanan .seealso `DMSlicedGetMatrix()`, `DMDASetBlockFills()` 18647c6ae99SBarry Smith @*/ 187*d71ae5a4SJacob Faibussowitsch PetscErrorCode DMSlicedSetBlockFills(DM dm, const PetscInt *dfill, const PetscInt *ofill) 188*d71ae5a4SJacob Faibussowitsch { 18947c6ae99SBarry Smith DM_Sliced *slice = (DM_Sliced *)dm->data; 19047c6ae99SBarry Smith 19147c6ae99SBarry Smith PetscFunctionBegin; 1920c010503SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1939566063dSJacob Faibussowitsch PetscCall(DMSlicedSetBlockFills_Private(slice->bs, dfill, &slice->dfill)); 1949566063dSJacob Faibussowitsch PetscCall(DMSlicedSetBlockFills_Private(slice->bs, ofill, &slice->ofill)); 19547c6ae99SBarry Smith PetscFunctionReturn(0); 19647c6ae99SBarry Smith } 19747c6ae99SBarry Smith 198*d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDestroy_Sliced(DM dm) 199*d71ae5a4SJacob Faibussowitsch { 20047c6ae99SBarry Smith DM_Sliced *slice = (DM_Sliced *)dm->data; 20147c6ae99SBarry Smith 20247c6ae99SBarry Smith PetscFunctionBegin; 2039566063dSJacob Faibussowitsch PetscCall(PetscFree(slice->ghosts)); 2049566063dSJacob Faibussowitsch if (slice->dfill) PetscCall(PetscFree3(slice->dfill, slice->dfill->i, slice->dfill->j)); 2059566063dSJacob Faibussowitsch if (slice->ofill) PetscCall(PetscFree3(slice->ofill, slice->ofill->i, slice->ofill->j)); 206435a35e8SMatthew G Knepley /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */ 2079566063dSJacob Faibussowitsch PetscCall(PetscFree(slice)); 20847c6ae99SBarry Smith PetscFunctionReturn(0); 20947c6ae99SBarry Smith } 21047c6ae99SBarry Smith 211*d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMCreateGlobalVector_Sliced(DM dm, Vec *gvec) 212*d71ae5a4SJacob Faibussowitsch { 21347c6ae99SBarry Smith DM_Sliced *slice = (DM_Sliced *)dm->data; 21447c6ae99SBarry Smith 21547c6ae99SBarry Smith PetscFunctionBegin; 21647c6ae99SBarry Smith PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 21747c6ae99SBarry Smith PetscValidPointer(gvec, 2); 218ea78f98cSLisandro Dalcin *gvec = NULL; 2199566063dSJacob Faibussowitsch PetscCall(VecCreateGhostBlock(PetscObjectComm((PetscObject)dm), slice->bs, slice->n * slice->bs, PETSC_DETERMINE, slice->Nghosts, slice->ghosts, gvec)); 2209566063dSJacob Faibussowitsch PetscCall(VecSetDM(*gvec, dm)); 22147c6ae99SBarry Smith PetscFunctionReturn(0); 22247c6ae99SBarry Smith } 22347c6ae99SBarry Smith 224*d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMGlobalToLocalBegin_Sliced(DM da, Vec g, InsertMode mode, Vec l) 225*d71ae5a4SJacob Faibussowitsch { 2263efe6655SBarry Smith PetscBool flg; 2273efe6655SBarry Smith 2283efe6655SBarry Smith PetscFunctionBegin; 2299566063dSJacob Faibussowitsch PetscCall(VecGhostIsLocalForm(g, l, &flg)); 23028b400f6SJacob Faibussowitsch PetscCheck(flg, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONG, "Local vector is not local form of global vector"); 2319566063dSJacob Faibussowitsch PetscCall(VecGhostUpdateEnd(g, mode, SCATTER_FORWARD)); 2329566063dSJacob Faibussowitsch PetscCall(VecGhostUpdateBegin(g, mode, SCATTER_FORWARD)); 2333efe6655SBarry Smith PetscFunctionReturn(0); 2343efe6655SBarry Smith } 2353efe6655SBarry Smith 236*d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMGlobalToLocalEnd_Sliced(DM da, Vec g, InsertMode mode, Vec l) 237*d71ae5a4SJacob Faibussowitsch { 2383efe6655SBarry Smith PetscBool flg; 2393efe6655SBarry Smith 2403efe6655SBarry Smith PetscFunctionBegin; 2419566063dSJacob Faibussowitsch PetscCall(VecGhostIsLocalForm(g, l, &flg)); 24228b400f6SJacob Faibussowitsch PetscCheck(flg, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONG, "Local vector is not local form of global vector"); 2439566063dSJacob Faibussowitsch PetscCall(VecGhostUpdateEnd(g, mode, SCATTER_FORWARD)); 2443efe6655SBarry Smith PetscFunctionReturn(0); 2453efe6655SBarry Smith } 2463efe6655SBarry Smith 2473e0de288SBarry Smith /*MC 2483e0de288SBarry Smith DMSLICED = "sliced" - A DM object that is used to manage data for a general graph. Uses VecCreateGhost() ghosted vectors for storing the fields 2493e0de288SBarry Smith 2503e0de288SBarry Smith See DMCreateSliced() for details. 2513e0de288SBarry Smith 2523e0de288SBarry Smith Level: intermediate 2533e0de288SBarry Smith 254db781477SPatrick Sanan .seealso: `DMType`, `DMCOMPOSITE`, `DMCreateSliced()`, `DMCreate()` 2553e0de288SBarry Smith M*/ 2563e0de288SBarry Smith 257*d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMCreate_Sliced(DM p) 258*d71ae5a4SJacob Faibussowitsch { 259a4121054SBarry Smith DM_Sliced *slice; 260a4121054SBarry Smith 261a4121054SBarry Smith PetscFunctionBegin; 2624dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&slice)); 263a4121054SBarry Smith p->data = slice; 264a4121054SBarry Smith 265a4121054SBarry Smith p->ops->createglobalvector = DMCreateGlobalVector_Sliced; 26625296bd5SBarry Smith p->ops->creatematrix = DMCreateMatrix_Sliced; 2673efe6655SBarry Smith p->ops->globaltolocalbegin = DMGlobalToLocalBegin_Sliced; 2683efe6655SBarry Smith p->ops->globaltolocalend = DMGlobalToLocalEnd_Sliced; 269a4121054SBarry Smith p->ops->destroy = DMDestroy_Sliced; 270a4121054SBarry Smith PetscFunctionReturn(0); 271a4121054SBarry Smith } 272a4121054SBarry Smith 27347c6ae99SBarry Smith /*@C 2740c010503SBarry Smith DMSlicedCreate - Creates a DM object, used to manage data for a unstructured problem 27547c6ae99SBarry Smith 276d083f849SBarry Smith Collective 2770c010503SBarry Smith 278d8d19677SJose E. Roman Input Parameters: 2793efe6655SBarry Smith + comm - the processors that will share the global vector 2803efe6655SBarry Smith . bs - the block size 2813efe6655SBarry Smith . nlocal - number of vector entries on this process 2823efe6655SBarry Smith . Nghosts - number of ghost points needed on this process 2833efe6655SBarry Smith . ghosts - global indices of all ghost points for this process 2843efe6655SBarry Smith . d_nnz - matrix preallocation information representing coupling within this process 2853efe6655SBarry Smith - o_nnz - matrix preallocation information representing coupling between this process and other processes 2860c010503SBarry Smith 2870c010503SBarry Smith Output Parameters: 2880c010503SBarry Smith . slice - the slice object 2890c010503SBarry Smith 2903efe6655SBarry Smith Notes: 2913efe6655SBarry Smith This DM does not support DMCreateLocalVector(), DMGlobalToLocalBegin(), and DMGlobalToLocalEnd() instead one directly uses 2923efe6655SBarry Smith VecGhostGetLocalForm() and VecGhostRestoreLocalForm() to access the local representation and VecGhostUpdateBegin() and VecGhostUpdateEnd() to update 2933efe6655SBarry Smith the ghost points. 2943efe6655SBarry Smith 2953efe6655SBarry Smith One can use DMGlobalToLocalBegin(), and DMGlobalToLocalEnd() instead of VecGhostUpdateBegin() and VecGhostUpdateEnd(). 2963efe6655SBarry Smith 2970c010503SBarry Smith Level: advanced 2980c010503SBarry Smith 299db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCreateGlobalVector()`, `DMSetType()`, `DMSLICED`, `DMSlicedSetGhosts()`, `DMSlicedSetPreallocation()`, `VecGhostUpdateBegin()`, `VecGhostUpdateEnd()`, 300db781477SPatrick Sanan `VecGhostGetLocalForm()`, `VecGhostRestoreLocalForm()` 3010c010503SBarry Smith 3020c010503SBarry Smith @*/ 303*d71ae5a4SJacob 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) 304*d71ae5a4SJacob Faibussowitsch { 3050c010503SBarry Smith PetscFunctionBegin; 306064a246eSJacob Faibussowitsch PetscValidPointer(dm, 8); 3079566063dSJacob Faibussowitsch PetscCall(DMCreate(comm, dm)); 3089566063dSJacob Faibussowitsch PetscCall(DMSetType(*dm, DMSLICED)); 3099566063dSJacob Faibussowitsch PetscCall(DMSlicedSetGhosts(*dm, bs, nlocal, Nghosts, ghosts)); 3101baa6e33SBarry Smith if (d_nnz) PetscCall(DMSlicedSetPreallocation(*dm, 0, d_nnz, 0, o_nnz)); 3110c010503SBarry Smith PetscFunctionReturn(0); 3120c010503SBarry Smith } 313