xref: /petsc/src/dm/impls/sliced/sliced.c (revision 28b400f66ebc7ae0049166a2294dfcd3df27e64b)
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 
16b412c318SBarry Smith PetscErrorCode  DMCreateMatrix_Sliced(DM dm, Mat *J)
1747c6ae99SBarry Smith {
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;
255f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PetscObjectComm((PetscObject)dm),J));
265f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(*J,slice->n*bs,slice->n*bs,PETSC_DETERMINE,PETSC_DETERMINE));
275f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetBlockSize(*J,bs));
285f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(*J,dm->mattype));
295f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSeqBAIJSetPreallocation(*J,bs,slice->d_nz,slice->d_nnz));
305f80ce2aSJacob Faibussowitsch   CHKERRQ(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. */
335f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectQueryFunction((PetscObject)*J,"MatMPIAIJSetPreallocation_C",&aij));
3447c6ae99SBarry Smith   if (!aij) {
355f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscObjectQueryFunction((PetscObject)*J,"MatSeqAIJSetPreallocation_C",&aij));
3647c6ae99SBarry Smith   }
3747c6ae99SBarry Smith   if (aij) {
3847c6ae99SBarry Smith     if (bs == 1) {
395f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSeqAIJSetPreallocation(*J,slice->d_nz,slice->d_nnz));
405f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMPIAIJSetPreallocation(*J,slice->d_nz,slice->d_nnz,slice->o_nz,slice->o_nnz));
4147c6ae99SBarry Smith     } else if (!slice->d_nnz) {
425f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSeqAIJSetPreallocation(*J,slice->d_nz*bs,NULL));
435f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMPIAIJSetPreallocation(*J,slice->d_nz*bs,NULL,slice->o_nz*bs,NULL));
4447c6ae99SBarry Smith     } else {
450c010503SBarry Smith       /* The user has provided preallocation per block-row, convert it to per scalar-row respecting DMSlicedSetBlockFills() if applicable */
465f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscMalloc2(slice->n*bs,&sd_nnz,(!!slice->o_nnz)*slice->n*bs,&so_nnz));
4747c6ae99SBarry Smith       for (i=0; i<slice->n*bs; i++) {
4847c6ae99SBarry Smith         sd_nnz[i] = (slice->d_nnz[i/bs]-1) * (slice->ofill ? slice->ofill->i[i%bs+1]-slice->ofill->i[i%bs] : bs)
4947c6ae99SBarry Smith                                            + (slice->dfill ? slice->dfill->i[i%bs+1]-slice->dfill->i[i%bs] : bs);
5047c6ae99SBarry Smith         if (so_nnz) {
5147c6ae99SBarry Smith           so_nnz[i] = slice->o_nnz[i/bs] * (slice->ofill ? slice->ofill->i[i%bs+1]-slice->ofill->i[i%bs] : bs);
5247c6ae99SBarry Smith         }
5347c6ae99SBarry Smith       }
545f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSeqAIJSetPreallocation(*J,slice->d_nz*bs,sd_nnz));
555f80ce2aSJacob Faibussowitsch       CHKERRQ(MatMPIAIJSetPreallocation(*J,slice->d_nz*bs,sd_nnz,slice->o_nz*bs,so_nnz));
565f80ce2aSJacob Faibussowitsch       CHKERRQ(PetscFree2(sd_nnz,so_nnz));
5747c6ae99SBarry Smith     }
5847c6ae99SBarry Smith   }
5947c6ae99SBarry Smith 
6047c6ae99SBarry Smith   /* Set up the local to global map.  For the scalar map, we have to translate to entry-wise indexing instead of block-wise. */
615f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(slice->n+slice->Nghosts,&globals));
625f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetOwnershipRange(*J,&rstart,NULL));
6345b6f7e9SBarry Smith   for (i=0; i<slice->n; i++) globals[i] = rstart/bs + i;
648865f1eaSKarl Rupp 
6545b6f7e9SBarry Smith   for (i=0; i<slice->Nghosts; i++) {
6645b6f7e9SBarry Smith     globals[slice->n+i] = slice->ghosts[i];
6747c6ae99SBarry Smith   }
685f80ce2aSJacob Faibussowitsch   CHKERRQ(ISLocalToGlobalMappingCreate(PETSC_COMM_SELF,bs,slice->n+slice->Nghosts,globals,PETSC_OWN_POINTER,&lmap));
695f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetLocalToGlobalMapping(*J,lmap,lmap));
705f80ce2aSJacob Faibussowitsch   CHKERRQ(ISLocalToGlobalMappingDestroy(&lmap));
715f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetDM(*J,dm));
7247c6ae99SBarry Smith   PetscFunctionReturn(0);
7347c6ae99SBarry Smith }
7447c6ae99SBarry Smith 
7547c6ae99SBarry Smith /*@C
760c010503SBarry Smith     DMSlicedSetGhosts - Sets the global indices of other processes elements that will
7747c6ae99SBarry Smith       be ghosts on this process
7847c6ae99SBarry Smith 
7947c6ae99SBarry Smith     Not Collective
8047c6ae99SBarry Smith 
8147c6ae99SBarry Smith     Input Parameters:
820c010503SBarry Smith +    slice - the DM object
8347c6ae99SBarry Smith .    bs - block size
8447c6ae99SBarry Smith .    nlocal - number of local (owned, non-ghost) blocks
8547c6ae99SBarry Smith .    Nghosts - number of ghost blocks on this process
8647c6ae99SBarry Smith -    ghosts - global indices of each ghost block
8747c6ae99SBarry Smith 
8847c6ae99SBarry Smith     Level: advanced
8947c6ae99SBarry Smith 
90240edb01SJed Brown .seealso DMDestroy(), DMCreateGlobalVector()
9147c6ae99SBarry Smith 
9247c6ae99SBarry Smith @*/
937087cfbeSBarry Smith PetscErrorCode  DMSlicedSetGhosts(DM dm,PetscInt bs,PetscInt nlocal,PetscInt Nghosts,const PetscInt ghosts[])
9447c6ae99SBarry Smith {
9547c6ae99SBarry Smith   DM_Sliced      *slice = (DM_Sliced*)dm->data;
9647c6ae99SBarry Smith 
9747c6ae99SBarry Smith   PetscFunctionBegin;
980c010503SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
995f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(slice->ghosts));
1005f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(Nghosts,&slice->ghosts));
1015f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscArraycpy(slice->ghosts,ghosts,Nghosts));
10247c6ae99SBarry Smith   slice->bs      = bs;
10347c6ae99SBarry Smith   slice->n       = nlocal;
10447c6ae99SBarry Smith   slice->Nghosts = Nghosts;
10547c6ae99SBarry Smith   PetscFunctionReturn(0);
10647c6ae99SBarry Smith }
10747c6ae99SBarry Smith 
10847c6ae99SBarry Smith /*@C
1090c010503SBarry Smith     DMSlicedSetPreallocation - sets the matrix memory preallocation for matrices computed by DMSliced
11047c6ae99SBarry Smith 
11147c6ae99SBarry Smith     Not Collective
11247c6ae99SBarry Smith 
11347c6ae99SBarry Smith     Input Parameters:
1140c010503SBarry Smith +    slice - the DM object
11547c6ae99SBarry Smith .    d_nz  - number of block nonzeros per block row in diagonal portion of local
11647c6ae99SBarry Smith            submatrix  (same for all local rows)
11747c6ae99SBarry Smith .    d_nnz - array containing the number of block nonzeros in the various block rows
11847c6ae99SBarry Smith            of the in diagonal portion of the local (possibly different for each block
119f77a5242SKarl Rupp            row) or NULL.
12047c6ae99SBarry Smith .    o_nz  - number of block nonzeros per block row in the off-diagonal portion of local
12147c6ae99SBarry Smith            submatrix (same for all local rows).
12247c6ae99SBarry Smith -    o_nnz - array containing the number of nonzeros in the various block rows of the
12347c6ae99SBarry Smith            off-diagonal portion of the local submatrix (possibly different for
124f77a5242SKarl Rupp            each block row) or NULL.
12547c6ae99SBarry Smith 
12647c6ae99SBarry Smith     Notes:
12747c6ae99SBarry Smith     See MatMPIBAIJSetPreallocation() for more details on preallocation.  If a scalar matrix (AIJ) is
1280c010503SBarry Smith     obtained with DMSlicedGetMatrix(), the correct preallocation will be set, respecting DMSlicedSetBlockFills().
12947c6ae99SBarry Smith 
13047c6ae99SBarry Smith     Level: advanced
13147c6ae99SBarry Smith 
132240edb01SJed Brown .seealso DMDestroy(), DMCreateGlobalVector(), MatMPIAIJSetPreallocation(),
1330c010503SBarry Smith          MatMPIBAIJSetPreallocation(), DMSlicedGetMatrix(), DMSlicedSetBlockFills()
13447c6ae99SBarry Smith 
13547c6ae99SBarry Smith @*/
1367087cfbeSBarry Smith PetscErrorCode  DMSlicedSetPreallocation(DM dm,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
13747c6ae99SBarry Smith {
13847c6ae99SBarry Smith   DM_Sliced *slice = (DM_Sliced*)dm->data;
13947c6ae99SBarry Smith 
14047c6ae99SBarry Smith   PetscFunctionBegin;
14147c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
14247c6ae99SBarry Smith   slice->d_nz  = d_nz;
14347c6ae99SBarry Smith   slice->d_nnz = (PetscInt*)d_nnz;
14447c6ae99SBarry Smith   slice->o_nz  = o_nz;
14547c6ae99SBarry Smith   slice->o_nnz = (PetscInt*)o_nnz;
14647c6ae99SBarry Smith   PetscFunctionReturn(0);
14747c6ae99SBarry Smith }
14847c6ae99SBarry Smith 
1490c010503SBarry Smith static PetscErrorCode DMSlicedSetBlockFills_Private(PetscInt bs,const PetscInt *fill,DMSlicedBlockFills **inf)
15047c6ae99SBarry Smith {
15147c6ae99SBarry Smith   PetscInt           i,j,nz,*fi,*fj;
1520c010503SBarry Smith   DMSlicedBlockFills *f;
15347c6ae99SBarry Smith 
15447c6ae99SBarry Smith   PetscFunctionBegin;
15547c6ae99SBarry Smith   PetscValidPointer(inf,3);
1565f80ce2aSJacob Faibussowitsch   if (*inf) CHKERRQ(PetscFree3(*inf,(*inf)->i,(*inf)->j));
15747c6ae99SBarry Smith   if (!fill) PetscFunctionReturn(0);
15847c6ae99SBarry Smith   for (i=0,nz=0; i<bs*bs; i++) if (fill[i]) nz++;
1595f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc3(1,&f,bs+1,&fi,nz,&fj));
16047c6ae99SBarry Smith   f->bs = bs;
16147c6ae99SBarry Smith   f->nz = nz;
16247c6ae99SBarry Smith   f->i  = fi;
16347c6ae99SBarry Smith   f->j  = fj;
16447c6ae99SBarry Smith   for (i=0,nz=0; i<bs; i++) {
16547c6ae99SBarry Smith     fi[i] = nz;
16647c6ae99SBarry Smith     for (j=0; j<bs; j++) if (fill[i*bs+j]) fj[nz++] = j;
16747c6ae99SBarry Smith   }
16847c6ae99SBarry Smith   fi[i] = nz;
16947c6ae99SBarry Smith   *inf  = f;
17047c6ae99SBarry Smith   PetscFunctionReturn(0);
17147c6ae99SBarry Smith }
17247c6ae99SBarry Smith 
17347c6ae99SBarry Smith /*@C
1740c010503SBarry Smith     DMSlicedSetBlockFills - Sets the fill pattern in each block for a multi-component problem
1750c010503SBarry Smith     of the matrix returned by DMSlicedGetMatrix().
17647c6ae99SBarry Smith 
177d083f849SBarry Smith     Logically Collective on dm
17847c6ae99SBarry Smith 
179d8d19677SJose E. Roman     Input Parameters:
1800c010503SBarry Smith +   sliced - the DM object
181f77a5242SKarl Rupp .   dfill - the fill pattern in the diagonal block (may be NULL, means use dense block)
18247c6ae99SBarry Smith -   ofill - the fill pattern in the off-diagonal blocks
18347c6ae99SBarry Smith 
18447c6ae99SBarry Smith     Notes:
18547c6ae99SBarry Smith     This only makes sense for multicomponent problems using scalar matrix formats (AIJ).
186aa219208SBarry Smith     See DMDASetBlockFills() for example usage.
18747c6ae99SBarry Smith 
18847c6ae99SBarry Smith     Level: advanced
18947c6ae99SBarry Smith 
190aa219208SBarry Smith .seealso DMSlicedGetMatrix(), DMDASetBlockFills()
19147c6ae99SBarry Smith @*/
1927087cfbeSBarry Smith PetscErrorCode  DMSlicedSetBlockFills(DM dm,const PetscInt *dfill,const PetscInt *ofill)
19347c6ae99SBarry Smith {
19447c6ae99SBarry Smith   DM_Sliced      *slice = (DM_Sliced*)dm->data;
19547c6ae99SBarry Smith 
19647c6ae99SBarry Smith   PetscFunctionBegin;
1970c010503SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1985f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSlicedSetBlockFills_Private(slice->bs,dfill,&slice->dfill));
1995f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSlicedSetBlockFills_Private(slice->bs,ofill,&slice->ofill));
20047c6ae99SBarry Smith   PetscFunctionReturn(0);
20147c6ae99SBarry Smith }
20247c6ae99SBarry Smith 
2033efe6655SBarry Smith static PetscErrorCode  DMDestroy_Sliced(DM dm)
20447c6ae99SBarry Smith {
20547c6ae99SBarry Smith   DM_Sliced      *slice = (DM_Sliced*)dm->data;
20647c6ae99SBarry Smith 
20747c6ae99SBarry Smith   PetscFunctionBegin;
2085f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(slice->ghosts));
2095f80ce2aSJacob Faibussowitsch   if (slice->dfill) CHKERRQ(PetscFree3(slice->dfill,slice->dfill->i,slice->dfill->j));
2105f80ce2aSJacob Faibussowitsch   if (slice->ofill) CHKERRQ(PetscFree3(slice->ofill,slice->ofill->i,slice->ofill->j));
211435a35e8SMatthew G Knepley   /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */
2125f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(slice));
21347c6ae99SBarry Smith   PetscFunctionReturn(0);
21447c6ae99SBarry Smith }
21547c6ae99SBarry Smith 
2163efe6655SBarry Smith static PetscErrorCode  DMCreateGlobalVector_Sliced(DM dm,Vec *gvec)
21747c6ae99SBarry Smith {
21847c6ae99SBarry Smith   DM_Sliced      *slice = (DM_Sliced*)dm->data;
21947c6ae99SBarry Smith 
22047c6ae99SBarry Smith   PetscFunctionBegin;
22147c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
22247c6ae99SBarry Smith   PetscValidPointer(gvec,2);
223ea78f98cSLisandro Dalcin   *gvec = NULL;
2245f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateGhostBlock(PetscObjectComm((PetscObject)dm),slice->bs,slice->n*slice->bs,PETSC_DETERMINE,slice->Nghosts,slice->ghosts,gvec));
2255f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetDM(*gvec,dm));
22647c6ae99SBarry Smith   PetscFunctionReturn(0);
22747c6ae99SBarry Smith }
22847c6ae99SBarry Smith 
2293efe6655SBarry Smith static PetscErrorCode  DMGlobalToLocalBegin_Sliced(DM da,Vec g,InsertMode mode,Vec l)
2303efe6655SBarry Smith {
2313efe6655SBarry Smith   PetscBool      flg;
2323efe6655SBarry Smith 
2333efe6655SBarry Smith   PetscFunctionBegin;
2345f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGhostIsLocalForm(g,l,&flg));
235*28b400f6SJacob Faibussowitsch   PetscCheck(flg,PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONG,"Local vector is not local form of global vector");
2365f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGhostUpdateEnd(g,mode,SCATTER_FORWARD));
2375f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGhostUpdateBegin(g,mode,SCATTER_FORWARD));
2383efe6655SBarry Smith   PetscFunctionReturn(0);
2393efe6655SBarry Smith }
2403efe6655SBarry Smith 
2413efe6655SBarry Smith static PetscErrorCode  DMGlobalToLocalEnd_Sliced(DM da,Vec g,InsertMode mode,Vec l)
2423efe6655SBarry Smith {
2433efe6655SBarry Smith   PetscBool      flg;
2443efe6655SBarry Smith 
2453efe6655SBarry Smith   PetscFunctionBegin;
2465f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGhostIsLocalForm(g,l,&flg));
247*28b400f6SJacob Faibussowitsch   PetscCheck(flg,PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONG,"Local vector is not local form of global vector");
2485f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGhostUpdateEnd(g,mode,SCATTER_FORWARD));
2493efe6655SBarry Smith   PetscFunctionReturn(0);
2503efe6655SBarry Smith }
2513efe6655SBarry Smith 
2523e0de288SBarry Smith /*MC
2533e0de288SBarry Smith    DMSLICED = "sliced" - A DM object that is used to manage data for a general graph. Uses VecCreateGhost() ghosted vectors for storing the fields
2543e0de288SBarry Smith 
2553e0de288SBarry Smith    See DMCreateSliced() for details.
2563e0de288SBarry Smith 
2573e0de288SBarry Smith   Level: intermediate
2583e0de288SBarry Smith 
2593e0de288SBarry Smith .seealso: DMType, DMCOMPOSITE, DMCreateSliced(), DMCreate()
2603e0de288SBarry Smith M*/
2613e0de288SBarry Smith 
2628cc058d9SJed Brown PETSC_EXTERN PetscErrorCode DMCreate_Sliced(DM p)
263a4121054SBarry Smith {
264a4121054SBarry Smith   DM_Sliced      *slice;
265a4121054SBarry Smith 
266a4121054SBarry Smith   PetscFunctionBegin;
2675f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscNewLog(p,&slice));
268a4121054SBarry Smith   p->data = slice;
269a4121054SBarry Smith 
270a4121054SBarry Smith   p->ops->createglobalvector = DMCreateGlobalVector_Sliced;
27125296bd5SBarry Smith   p->ops->creatematrix       = DMCreateMatrix_Sliced;
2723efe6655SBarry Smith   p->ops->globaltolocalbegin = DMGlobalToLocalBegin_Sliced;
2733efe6655SBarry Smith   p->ops->globaltolocalend   = DMGlobalToLocalEnd_Sliced;
274a4121054SBarry Smith   p->ops->destroy            = DMDestroy_Sliced;
275a4121054SBarry Smith   PetscFunctionReturn(0);
276a4121054SBarry Smith }
277a4121054SBarry Smith 
27847c6ae99SBarry Smith /*@C
2790c010503SBarry Smith     DMSlicedCreate - Creates a DM object, used to manage data for a unstructured problem
28047c6ae99SBarry Smith 
281d083f849SBarry Smith     Collective
2820c010503SBarry Smith 
283d8d19677SJose E. Roman     Input Parameters:
2843efe6655SBarry Smith +   comm - the processors that will share the global vector
2853efe6655SBarry Smith .   bs - the block size
2863efe6655SBarry Smith .   nlocal - number of vector entries on this process
2873efe6655SBarry Smith .   Nghosts - number of ghost points needed on this process
2883efe6655SBarry Smith .   ghosts - global indices of all ghost points for this process
2893efe6655SBarry Smith .   d_nnz - matrix preallocation information representing coupling within this process
2903efe6655SBarry Smith -   o_nnz - matrix preallocation information representing coupling between this process and other processes
2910c010503SBarry Smith 
2920c010503SBarry Smith     Output Parameters:
2930c010503SBarry Smith .   slice - the slice object
2940c010503SBarry Smith 
2953efe6655SBarry Smith     Notes:
2963efe6655SBarry Smith         This DM does not support DMCreateLocalVector(), DMGlobalToLocalBegin(), and DMGlobalToLocalEnd() instead one directly uses
2973efe6655SBarry Smith         VecGhostGetLocalForm() and VecGhostRestoreLocalForm() to access the local representation and VecGhostUpdateBegin() and VecGhostUpdateEnd() to update
2983efe6655SBarry Smith         the ghost points.
2993efe6655SBarry Smith 
3003efe6655SBarry Smith         One can use DMGlobalToLocalBegin(), and DMGlobalToLocalEnd() instead of VecGhostUpdateBegin() and VecGhostUpdateEnd().
3013efe6655SBarry Smith 
3020c010503SBarry Smith     Level: advanced
3030c010503SBarry Smith 
3043e0de288SBarry Smith .seealso DMDestroy(), DMCreateGlobalVector(), DMSetType(), DMSLICED, DMSlicedSetGhosts(), DMSlicedSetPreallocation(), VecGhostUpdateBegin(), VecGhostUpdateEnd(),
3053e0de288SBarry Smith          VecGhostGetLocalForm(), VecGhostRestoreLocalForm()
3060c010503SBarry Smith 
3070c010503SBarry Smith @*/
3083efe6655SBarry Smith PetscErrorCode  DMSlicedCreate(MPI_Comm comm,PetscInt bs,PetscInt nlocal,PetscInt Nghosts,const PetscInt ghosts[], const PetscInt d_nnz[],const PetscInt o_nnz[],DM *dm)
3090c010503SBarry Smith {
3100c010503SBarry Smith   PetscFunctionBegin;
311064a246eSJacob Faibussowitsch   PetscValidPointer(dm,8);
3125f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreate(comm,dm));
3135f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetType(*dm,DMSLICED));
3145f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSlicedSetGhosts(*dm,bs,nlocal,Nghosts,ghosts));
3153efe6655SBarry Smith   if (d_nnz) {
3165f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSlicedSetPreallocation(*dm,0, d_nnz,0,o_nnz));
3173efe6655SBarry Smith   }
3180c010503SBarry Smith   PetscFunctionReturn(0);
3190c010503SBarry Smith }
320