xref: /petsc/src/dm/impls/sliced/sliced.c (revision 4dfa11a44d5adf2389f1d3acbc8f3c1116dc6c3a)
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 
169371c9d4SSatish Balay PetscErrorCode DMCreateMatrix_Sliced(DM dm, Mat *J) {
1747c6ae99SBarry Smith   PetscInt              *globals, *sd_nnz, *so_nnz, rstart, bs, i;
1845b6f7e9SBarry Smith   ISLocalToGlobalMapping lmap;
19f77a5242SKarl Rupp   void (*aij)(void) = NULL;
2047c6ae99SBarry Smith   DM_Sliced *slice  = (DM_Sliced *)dm->data;
2147c6ae99SBarry Smith 
2247c6ae99SBarry Smith   PetscFunctionBegin;
2347c6ae99SBarry Smith   bs = slice->bs;
249566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)dm), J));
259566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(*J, slice->n * bs, slice->n * bs, PETSC_DETERMINE, PETSC_DETERMINE));
269566063dSJacob Faibussowitsch   PetscCall(MatSetBlockSize(*J, bs));
279566063dSJacob Faibussowitsch   PetscCall(MatSetType(*J, dm->mattype));
289566063dSJacob Faibussowitsch   PetscCall(MatSeqBAIJSetPreallocation(*J, bs, slice->d_nz, slice->d_nnz));
299566063dSJacob Faibussowitsch   PetscCall(MatMPIBAIJSetPreallocation(*J, bs, slice->d_nz, slice->d_nnz, slice->o_nz, slice->o_nnz));
3047c6ae99SBarry Smith   /* In general, we have to do extra work to preallocate for scalar (AIJ) matrices so we check whether it will do any
3147c6ae99SBarry Smith   * good before going on with it. */
329566063dSJacob Faibussowitsch   PetscCall(PetscObjectQueryFunction((PetscObject)*J, "MatMPIAIJSetPreallocation_C", &aij));
3348a46eb9SPierre Jolivet   if (!aij) PetscCall(PetscObjectQueryFunction((PetscObject)*J, "MatSeqAIJSetPreallocation_C", &aij));
3447c6ae99SBarry Smith   if (aij) {
3547c6ae99SBarry Smith     if (bs == 1) {
369566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJSetPreallocation(*J, slice->d_nz, slice->d_nnz));
379566063dSJacob Faibussowitsch       PetscCall(MatMPIAIJSetPreallocation(*J, slice->d_nz, slice->d_nnz, slice->o_nz, slice->o_nnz));
3847c6ae99SBarry Smith     } else if (!slice->d_nnz) {
399566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJSetPreallocation(*J, slice->d_nz * bs, NULL));
409566063dSJacob Faibussowitsch       PetscCall(MatMPIAIJSetPreallocation(*J, slice->d_nz * bs, NULL, slice->o_nz * bs, NULL));
4147c6ae99SBarry Smith     } else {
420c010503SBarry Smith       /* The user has provided preallocation per block-row, convert it to per scalar-row respecting DMSlicedSetBlockFills() if applicable */
439566063dSJacob Faibussowitsch       PetscCall(PetscMalloc2(slice->n * bs, &sd_nnz, (!!slice->o_nnz) * slice->n * bs, &so_nnz));
4447c6ae99SBarry Smith       for (i = 0; i < slice->n * bs; i++) {
459371c9d4SSatish 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);
46ad540459SPierre 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);
4747c6ae99SBarry Smith       }
489566063dSJacob Faibussowitsch       PetscCall(MatSeqAIJSetPreallocation(*J, slice->d_nz * bs, sd_nnz));
499566063dSJacob Faibussowitsch       PetscCall(MatMPIAIJSetPreallocation(*J, slice->d_nz * bs, sd_nnz, slice->o_nz * bs, so_nnz));
509566063dSJacob Faibussowitsch       PetscCall(PetscFree2(sd_nnz, so_nnz));
5147c6ae99SBarry Smith     }
5247c6ae99SBarry Smith   }
5347c6ae99SBarry Smith 
5447c6ae99SBarry Smith   /* Set up the local to global map.  For the scalar map, we have to translate to entry-wise indexing instead of block-wise. */
559566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(slice->n + slice->Nghosts, &globals));
569566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(*J, &rstart, NULL));
5745b6f7e9SBarry Smith   for (i = 0; i < slice->n; i++) globals[i] = rstart / bs + i;
588865f1eaSKarl Rupp 
59ad540459SPierre Jolivet   for (i = 0; i < slice->Nghosts; i++) globals[slice->n + i] = slice->ghosts[i];
609566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingCreate(PETSC_COMM_SELF, bs, slice->n + slice->Nghosts, globals, PETSC_OWN_POINTER, &lmap));
619566063dSJacob Faibussowitsch   PetscCall(MatSetLocalToGlobalMapping(*J, lmap, lmap));
629566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingDestroy(&lmap));
639566063dSJacob Faibussowitsch   PetscCall(MatSetDM(*J, dm));
6447c6ae99SBarry Smith   PetscFunctionReturn(0);
6547c6ae99SBarry Smith }
6647c6ae99SBarry Smith 
6747c6ae99SBarry Smith /*@C
680c010503SBarry Smith     DMSlicedSetGhosts - Sets the global indices of other processes elements that will
6947c6ae99SBarry Smith       be ghosts on this process
7047c6ae99SBarry Smith 
7147c6ae99SBarry Smith     Not Collective
7247c6ae99SBarry Smith 
7347c6ae99SBarry Smith     Input Parameters:
740c010503SBarry Smith +    slice - the DM object
7547c6ae99SBarry Smith .    bs - block size
7647c6ae99SBarry Smith .    nlocal - number of local (owned, non-ghost) blocks
7747c6ae99SBarry Smith .    Nghosts - number of ghost blocks on this process
7847c6ae99SBarry Smith -    ghosts - global indices of each ghost block
7947c6ae99SBarry Smith 
8047c6ae99SBarry Smith     Level: advanced
8147c6ae99SBarry Smith 
82db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCreateGlobalVector()`
8347c6ae99SBarry Smith 
8447c6ae99SBarry Smith @*/
859371c9d4SSatish Balay PetscErrorCode DMSlicedSetGhosts(DM dm, PetscInt bs, PetscInt nlocal, PetscInt Nghosts, const PetscInt ghosts[]) {
8647c6ae99SBarry Smith   DM_Sliced *slice = (DM_Sliced *)dm->data;
8747c6ae99SBarry Smith 
8847c6ae99SBarry Smith   PetscFunctionBegin;
890c010503SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
909566063dSJacob Faibussowitsch   PetscCall(PetscFree(slice->ghosts));
919566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(Nghosts, &slice->ghosts));
929566063dSJacob Faibussowitsch   PetscCall(PetscArraycpy(slice->ghosts, ghosts, Nghosts));
9347c6ae99SBarry Smith   slice->bs      = bs;
9447c6ae99SBarry Smith   slice->n       = nlocal;
9547c6ae99SBarry Smith   slice->Nghosts = Nghosts;
9647c6ae99SBarry Smith   PetscFunctionReturn(0);
9747c6ae99SBarry Smith }
9847c6ae99SBarry Smith 
9947c6ae99SBarry Smith /*@C
1000c010503SBarry Smith     DMSlicedSetPreallocation - sets the matrix memory preallocation for matrices computed by DMSliced
10147c6ae99SBarry Smith 
10247c6ae99SBarry Smith     Not Collective
10347c6ae99SBarry Smith 
10447c6ae99SBarry Smith     Input Parameters:
1050c010503SBarry Smith +    slice - the DM object
10647c6ae99SBarry Smith .    d_nz  - number of block nonzeros per block row in diagonal portion of local
10747c6ae99SBarry Smith            submatrix  (same for all local rows)
10847c6ae99SBarry Smith .    d_nnz - array containing the number of block nonzeros in the various block rows
10947c6ae99SBarry Smith            of the in diagonal portion of the local (possibly different for each block
110f77a5242SKarl Rupp            row) or NULL.
11147c6ae99SBarry Smith .    o_nz  - number of block nonzeros per block row in the off-diagonal portion of local
11247c6ae99SBarry Smith            submatrix (same for all local rows).
11347c6ae99SBarry Smith -    o_nnz - array containing the number of nonzeros in the various block rows of the
11447c6ae99SBarry Smith            off-diagonal portion of the local submatrix (possibly different for
115f77a5242SKarl Rupp            each block row) or NULL.
11647c6ae99SBarry Smith 
11747c6ae99SBarry Smith     Notes:
11847c6ae99SBarry Smith     See MatMPIBAIJSetPreallocation() for more details on preallocation.  If a scalar matrix (AIJ) is
1190c010503SBarry Smith     obtained with DMSlicedGetMatrix(), the correct preallocation will be set, respecting DMSlicedSetBlockFills().
12047c6ae99SBarry Smith 
12147c6ae99SBarry Smith     Level: advanced
12247c6ae99SBarry Smith 
123db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCreateGlobalVector()`, `MatMPIAIJSetPreallocation()`,
124db781477SPatrick Sanan          `MatMPIBAIJSetPreallocation()`, `DMSlicedGetMatrix()`, `DMSlicedSetBlockFills()`
12547c6ae99SBarry Smith 
12647c6ae99SBarry Smith @*/
1279371c9d4SSatish Balay PetscErrorCode DMSlicedSetPreallocation(DM dm, PetscInt d_nz, const PetscInt d_nnz[], PetscInt o_nz, const PetscInt o_nnz[]) {
12847c6ae99SBarry Smith   DM_Sliced *slice = (DM_Sliced *)dm->data;
12947c6ae99SBarry Smith 
13047c6ae99SBarry Smith   PetscFunctionBegin;
13147c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
13247c6ae99SBarry Smith   slice->d_nz  = d_nz;
13347c6ae99SBarry Smith   slice->d_nnz = (PetscInt *)d_nnz;
13447c6ae99SBarry Smith   slice->o_nz  = o_nz;
13547c6ae99SBarry Smith   slice->o_nnz = (PetscInt *)o_nnz;
13647c6ae99SBarry Smith   PetscFunctionReturn(0);
13747c6ae99SBarry Smith }
13847c6ae99SBarry Smith 
1399371c9d4SSatish Balay static PetscErrorCode DMSlicedSetBlockFills_Private(PetscInt bs, const PetscInt *fill, DMSlicedBlockFills **inf) {
14047c6ae99SBarry Smith   PetscInt            i, j, nz, *fi, *fj;
1410c010503SBarry Smith   DMSlicedBlockFills *f;
14247c6ae99SBarry Smith 
14347c6ae99SBarry Smith   PetscFunctionBegin;
14447c6ae99SBarry Smith   PetscValidPointer(inf, 3);
1459566063dSJacob Faibussowitsch   if (*inf) PetscCall(PetscFree3(*inf, (*inf)->i, (*inf)->j));
14647c6ae99SBarry Smith   if (!fill) PetscFunctionReturn(0);
1479371c9d4SSatish Balay   for (i = 0, nz = 0; i < bs * bs; i++)
1489371c9d4SSatish Balay     if (fill[i]) nz++;
1499566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(1, &f, bs + 1, &fi, nz, &fj));
15047c6ae99SBarry Smith   f->bs = bs;
15147c6ae99SBarry Smith   f->nz = nz;
15247c6ae99SBarry Smith   f->i  = fi;
15347c6ae99SBarry Smith   f->j  = fj;
15447c6ae99SBarry Smith   for (i = 0, nz = 0; i < bs; i++) {
15547c6ae99SBarry Smith     fi[i] = nz;
1569371c9d4SSatish Balay     for (j = 0; j < bs; j++)
1579371c9d4SSatish Balay       if (fill[i * bs + j]) fj[nz++] = j;
15847c6ae99SBarry Smith   }
15947c6ae99SBarry Smith   fi[i] = nz;
16047c6ae99SBarry Smith   *inf  = f;
16147c6ae99SBarry Smith   PetscFunctionReturn(0);
16247c6ae99SBarry Smith }
16347c6ae99SBarry Smith 
16447c6ae99SBarry Smith /*@C
1650c010503SBarry Smith     DMSlicedSetBlockFills - Sets the fill pattern in each block for a multi-component problem
1660c010503SBarry Smith     of the matrix returned by DMSlicedGetMatrix().
16747c6ae99SBarry Smith 
168d083f849SBarry Smith     Logically Collective on dm
16947c6ae99SBarry Smith 
170d8d19677SJose E. Roman     Input Parameters:
1710c010503SBarry Smith +   sliced - the DM object
172f77a5242SKarl Rupp .   dfill - the fill pattern in the diagonal block (may be NULL, means use dense block)
17347c6ae99SBarry Smith -   ofill - the fill pattern in the off-diagonal blocks
17447c6ae99SBarry Smith 
17547c6ae99SBarry Smith     Notes:
17647c6ae99SBarry Smith     This only makes sense for multicomponent problems using scalar matrix formats (AIJ).
177aa219208SBarry Smith     See DMDASetBlockFills() for example usage.
17847c6ae99SBarry Smith 
17947c6ae99SBarry Smith     Level: advanced
18047c6ae99SBarry Smith 
181db781477SPatrick Sanan .seealso `DMSlicedGetMatrix()`, `DMDASetBlockFills()`
18247c6ae99SBarry Smith @*/
1839371c9d4SSatish Balay PetscErrorCode DMSlicedSetBlockFills(DM dm, const PetscInt *dfill, const PetscInt *ofill) {
18447c6ae99SBarry Smith   DM_Sliced *slice = (DM_Sliced *)dm->data;
18547c6ae99SBarry Smith 
18647c6ae99SBarry Smith   PetscFunctionBegin;
1870c010503SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1889566063dSJacob Faibussowitsch   PetscCall(DMSlicedSetBlockFills_Private(slice->bs, dfill, &slice->dfill));
1899566063dSJacob Faibussowitsch   PetscCall(DMSlicedSetBlockFills_Private(slice->bs, ofill, &slice->ofill));
19047c6ae99SBarry Smith   PetscFunctionReturn(0);
19147c6ae99SBarry Smith }
19247c6ae99SBarry Smith 
1939371c9d4SSatish Balay static PetscErrorCode DMDestroy_Sliced(DM dm) {
19447c6ae99SBarry Smith   DM_Sliced *slice = (DM_Sliced *)dm->data;
19547c6ae99SBarry Smith 
19647c6ae99SBarry Smith   PetscFunctionBegin;
1979566063dSJacob Faibussowitsch   PetscCall(PetscFree(slice->ghosts));
1989566063dSJacob Faibussowitsch   if (slice->dfill) PetscCall(PetscFree3(slice->dfill, slice->dfill->i, slice->dfill->j));
1999566063dSJacob Faibussowitsch   if (slice->ofill) PetscCall(PetscFree3(slice->ofill, slice->ofill->i, slice->ofill->j));
200435a35e8SMatthew G Knepley   /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */
2019566063dSJacob Faibussowitsch   PetscCall(PetscFree(slice));
20247c6ae99SBarry Smith   PetscFunctionReturn(0);
20347c6ae99SBarry Smith }
20447c6ae99SBarry Smith 
2059371c9d4SSatish Balay static PetscErrorCode DMCreateGlobalVector_Sliced(DM dm, Vec *gvec) {
20647c6ae99SBarry Smith   DM_Sliced *slice = (DM_Sliced *)dm->data;
20747c6ae99SBarry Smith 
20847c6ae99SBarry Smith   PetscFunctionBegin;
20947c6ae99SBarry Smith   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
21047c6ae99SBarry Smith   PetscValidPointer(gvec, 2);
211ea78f98cSLisandro Dalcin   *gvec = NULL;
2129566063dSJacob Faibussowitsch   PetscCall(VecCreateGhostBlock(PetscObjectComm((PetscObject)dm), slice->bs, slice->n * slice->bs, PETSC_DETERMINE, slice->Nghosts, slice->ghosts, gvec));
2139566063dSJacob Faibussowitsch   PetscCall(VecSetDM(*gvec, dm));
21447c6ae99SBarry Smith   PetscFunctionReturn(0);
21547c6ae99SBarry Smith }
21647c6ae99SBarry Smith 
2179371c9d4SSatish Balay static PetscErrorCode DMGlobalToLocalBegin_Sliced(DM da, Vec g, InsertMode mode, Vec l) {
2183efe6655SBarry Smith   PetscBool flg;
2193efe6655SBarry Smith 
2203efe6655SBarry Smith   PetscFunctionBegin;
2219566063dSJacob Faibussowitsch   PetscCall(VecGhostIsLocalForm(g, l, &flg));
22228b400f6SJacob Faibussowitsch   PetscCheck(flg, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONG, "Local vector is not local form of global vector");
2239566063dSJacob Faibussowitsch   PetscCall(VecGhostUpdateEnd(g, mode, SCATTER_FORWARD));
2249566063dSJacob Faibussowitsch   PetscCall(VecGhostUpdateBegin(g, mode, SCATTER_FORWARD));
2253efe6655SBarry Smith   PetscFunctionReturn(0);
2263efe6655SBarry Smith }
2273efe6655SBarry Smith 
2289371c9d4SSatish Balay static PetscErrorCode DMGlobalToLocalEnd_Sliced(DM da, Vec g, InsertMode mode, Vec l) {
2293efe6655SBarry Smith   PetscBool flg;
2303efe6655SBarry Smith 
2313efe6655SBarry Smith   PetscFunctionBegin;
2329566063dSJacob Faibussowitsch   PetscCall(VecGhostIsLocalForm(g, l, &flg));
23328b400f6SJacob Faibussowitsch   PetscCheck(flg, PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_WRONG, "Local vector is not local form of global vector");
2349566063dSJacob Faibussowitsch   PetscCall(VecGhostUpdateEnd(g, mode, SCATTER_FORWARD));
2353efe6655SBarry Smith   PetscFunctionReturn(0);
2363efe6655SBarry Smith }
2373efe6655SBarry Smith 
2383e0de288SBarry Smith /*MC
2393e0de288SBarry Smith    DMSLICED = "sliced" - A DM object that is used to manage data for a general graph. Uses VecCreateGhost() ghosted vectors for storing the fields
2403e0de288SBarry Smith 
2413e0de288SBarry Smith    See DMCreateSliced() for details.
2423e0de288SBarry Smith 
2433e0de288SBarry Smith   Level: intermediate
2443e0de288SBarry Smith 
245db781477SPatrick Sanan .seealso: `DMType`, `DMCOMPOSITE`, `DMCreateSliced()`, `DMCreate()`
2463e0de288SBarry Smith M*/
2473e0de288SBarry Smith 
2489371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode DMCreate_Sliced(DM p) {
249a4121054SBarry Smith   DM_Sliced *slice;
250a4121054SBarry Smith 
251a4121054SBarry Smith   PetscFunctionBegin;
252*4dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&slice));
253a4121054SBarry Smith   p->data = slice;
254a4121054SBarry Smith 
255a4121054SBarry Smith   p->ops->createglobalvector = DMCreateGlobalVector_Sliced;
25625296bd5SBarry Smith   p->ops->creatematrix       = DMCreateMatrix_Sliced;
2573efe6655SBarry Smith   p->ops->globaltolocalbegin = DMGlobalToLocalBegin_Sliced;
2583efe6655SBarry Smith   p->ops->globaltolocalend   = DMGlobalToLocalEnd_Sliced;
259a4121054SBarry Smith   p->ops->destroy            = DMDestroy_Sliced;
260a4121054SBarry Smith   PetscFunctionReturn(0);
261a4121054SBarry Smith }
262a4121054SBarry Smith 
26347c6ae99SBarry Smith /*@C
2640c010503SBarry Smith     DMSlicedCreate - Creates a DM object, used to manage data for a unstructured problem
26547c6ae99SBarry Smith 
266d083f849SBarry Smith     Collective
2670c010503SBarry Smith 
268d8d19677SJose E. Roman     Input Parameters:
2693efe6655SBarry Smith +   comm - the processors that will share the global vector
2703efe6655SBarry Smith .   bs - the block size
2713efe6655SBarry Smith .   nlocal - number of vector entries on this process
2723efe6655SBarry Smith .   Nghosts - number of ghost points needed on this process
2733efe6655SBarry Smith .   ghosts - global indices of all ghost points for this process
2743efe6655SBarry Smith .   d_nnz - matrix preallocation information representing coupling within this process
2753efe6655SBarry Smith -   o_nnz - matrix preallocation information representing coupling between this process and other processes
2760c010503SBarry Smith 
2770c010503SBarry Smith     Output Parameters:
2780c010503SBarry Smith .   slice - the slice object
2790c010503SBarry Smith 
2803efe6655SBarry Smith     Notes:
2813efe6655SBarry Smith         This DM does not support DMCreateLocalVector(), DMGlobalToLocalBegin(), and DMGlobalToLocalEnd() instead one directly uses
2823efe6655SBarry Smith         VecGhostGetLocalForm() and VecGhostRestoreLocalForm() to access the local representation and VecGhostUpdateBegin() and VecGhostUpdateEnd() to update
2833efe6655SBarry Smith         the ghost points.
2843efe6655SBarry Smith 
2853efe6655SBarry Smith         One can use DMGlobalToLocalBegin(), and DMGlobalToLocalEnd() instead of VecGhostUpdateBegin() and VecGhostUpdateEnd().
2863efe6655SBarry Smith 
2870c010503SBarry Smith     Level: advanced
2880c010503SBarry Smith 
289db781477SPatrick Sanan .seealso `DMDestroy()`, `DMCreateGlobalVector()`, `DMSetType()`, `DMSLICED`, `DMSlicedSetGhosts()`, `DMSlicedSetPreallocation()`, `VecGhostUpdateBegin()`, `VecGhostUpdateEnd()`,
290db781477SPatrick Sanan          `VecGhostGetLocalForm()`, `VecGhostRestoreLocalForm()`
2910c010503SBarry Smith 
2920c010503SBarry Smith @*/
2939371c9d4SSatish Balay PetscErrorCode DMSlicedCreate(MPI_Comm comm, PetscInt bs, PetscInt nlocal, PetscInt Nghosts, const PetscInt ghosts[], const PetscInt d_nnz[], const PetscInt o_nnz[], DM *dm) {
2940c010503SBarry Smith   PetscFunctionBegin;
295064a246eSJacob Faibussowitsch   PetscValidPointer(dm, 8);
2969566063dSJacob Faibussowitsch   PetscCall(DMCreate(comm, dm));
2979566063dSJacob Faibussowitsch   PetscCall(DMSetType(*dm, DMSLICED));
2989566063dSJacob Faibussowitsch   PetscCall(DMSlicedSetGhosts(*dm, bs, nlocal, Nghosts, ghosts));
2991baa6e33SBarry Smith   if (d_nnz) PetscCall(DMSlicedSetPreallocation(*dm, 0, d_nnz, 0, o_nnz));
3000c010503SBarry Smith   PetscFunctionReturn(0);
3010c010503SBarry Smith }
302