xref: /petsc/src/dm/impls/sliced/sliced.c (revision f73d5cc4fab9985badaa1a29b152f456bf59ad34)
1c6db04a5SJed Brown #include <petscdmsliced.h>      /*I      "petscdmsliced.h" I*/
2c6db04a5SJed Brown #include <petscmat.h>           /*I      "petscmat.h"      I*/
3b45d2f2cSJed Brown #include <petsc-private/dmimpl.h>     /*I      "petscdm.h"       I*/
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 
1647c6ae99SBarry Smith #undef __FUNCT__
17950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_Sliced"
18950540a4SJed Brown PetscErrorCode  DMCreateMatrix_Sliced(DM dm, const MatType mtype,Mat *J)
1947c6ae99SBarry Smith {
2047c6ae99SBarry Smith   PetscErrorCode         ierr;
2147c6ae99SBarry Smith   PetscInt               *globals,*sd_nnz,*so_nnz,rstart,bs,i;
2247c6ae99SBarry Smith   ISLocalToGlobalMapping lmap,blmap;
2347c6ae99SBarry Smith   void                   (*aij)(void) = PETSC_NULL;
2447c6ae99SBarry Smith   DM_Sliced              *slice = (DM_Sliced*)dm->data;
2547c6ae99SBarry Smith 
2647c6ae99SBarry Smith   PetscFunctionBegin;
2747c6ae99SBarry Smith   bs = slice->bs;
2847c6ae99SBarry Smith   ierr = MatCreate(((PetscObject)dm)->comm,J);CHKERRQ(ierr);
2947c6ae99SBarry Smith   ierr = MatSetSizes(*J,slice->n*bs,slice->n*bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
30*f73d5cc4SBarry Smith   ierr = MatSetBlockSize(*J,bs);CHKERRQ(ierr);
3147c6ae99SBarry Smith   ierr = MatSetType(*J,mtype);CHKERRQ(ierr);
3247c6ae99SBarry Smith   ierr = MatSeqBAIJSetPreallocation(*J,bs,slice->d_nz,slice->d_nnz);CHKERRQ(ierr);
3347c6ae99SBarry Smith   ierr = MatMPIBAIJSetPreallocation(*J,bs,slice->d_nz,slice->d_nnz,slice->o_nz,slice->o_nnz);CHKERRQ(ierr);
3447c6ae99SBarry Smith   /* In general, we have to do extra work to preallocate for scalar (AIJ) matrices so we check whether it will do any
3547c6ae99SBarry Smith   * good before going on with it. */
3647c6ae99SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)*J,"MatMPIAIJSetPreallocation_C",&aij);CHKERRQ(ierr);
3747c6ae99SBarry Smith   if (!aij) {
3847c6ae99SBarry Smith     ierr = PetscObjectQueryFunction((PetscObject)*J,"MatSeqAIJSetPreallocation_C",&aij);CHKERRQ(ierr);
3947c6ae99SBarry Smith   }
4047c6ae99SBarry Smith   if (aij) {
4147c6ae99SBarry Smith     if (bs == 1) {
4247c6ae99SBarry Smith       ierr = MatSeqAIJSetPreallocation(*J,slice->d_nz,slice->d_nnz);CHKERRQ(ierr);
4347c6ae99SBarry Smith       ierr = MatMPIAIJSetPreallocation(*J,slice->d_nz,slice->d_nnz,slice->o_nz,slice->o_nnz);CHKERRQ(ierr);
4447c6ae99SBarry Smith     } else if (!slice->d_nnz) {
4547c6ae99SBarry Smith       ierr = MatSeqAIJSetPreallocation(*J,slice->d_nz*bs,PETSC_NULL);CHKERRQ(ierr);
4647c6ae99SBarry Smith       ierr = MatMPIAIJSetPreallocation(*J,slice->d_nz*bs,PETSC_NULL,slice->o_nz*bs,PETSC_NULL);CHKERRQ(ierr);
4747c6ae99SBarry Smith     } else {
480c010503SBarry Smith       /* The user has provided preallocation per block-row, convert it to per scalar-row respecting DMSlicedSetBlockFills() if applicable */
4947c6ae99SBarry Smith       ierr = PetscMalloc2(slice->n*bs,PetscInt,&sd_nnz,(!!slice->o_nnz)*slice->n*bs,PetscInt,&so_nnz);CHKERRQ(ierr);
5047c6ae99SBarry Smith       for (i=0; i<slice->n*bs; i++) {
5147c6ae99SBarry Smith         sd_nnz[i] = (slice->d_nnz[i/bs]-1) * (slice->ofill ? slice->ofill->i[i%bs+1]-slice->ofill->i[i%bs] : bs)
5247c6ae99SBarry Smith                                            + (slice->dfill ? slice->dfill->i[i%bs+1]-slice->dfill->i[i%bs] : bs);
5347c6ae99SBarry Smith         if (so_nnz) {
5447c6ae99SBarry Smith           so_nnz[i] = slice->o_nnz[i/bs] * (slice->ofill ? slice->ofill->i[i%bs+1]-slice->ofill->i[i%bs] : bs);
5547c6ae99SBarry Smith         }
5647c6ae99SBarry Smith       }
5747c6ae99SBarry Smith       ierr = MatSeqAIJSetPreallocation(*J,slice->d_nz*bs,sd_nnz);CHKERRQ(ierr);
5847c6ae99SBarry Smith       ierr = MatMPIAIJSetPreallocation(*J,slice->d_nz*bs,sd_nnz,slice->o_nz*bs,so_nnz);CHKERRQ(ierr);
5947c6ae99SBarry Smith       ierr = PetscFree2(sd_nnz,so_nnz);CHKERRQ(ierr);
6047c6ae99SBarry Smith     }
6147c6ae99SBarry Smith   }
6247c6ae99SBarry Smith 
6347c6ae99SBarry Smith   /* Set up the local to global map.  For the scalar map, we have to translate to entry-wise indexing instead of block-wise. */
6447c6ae99SBarry Smith   ierr = PetscMalloc((slice->n+slice->Nghosts)*bs*sizeof(PetscInt),&globals);CHKERRQ(ierr);
6547c6ae99SBarry Smith   ierr = MatGetOwnershipRange(*J,&rstart,PETSC_NULL);CHKERRQ(ierr);
6647c6ae99SBarry Smith   for (i=0; i<slice->n*bs; i++) {
6747c6ae99SBarry Smith     globals[i] = rstart + i;
6847c6ae99SBarry Smith   }
6947c6ae99SBarry Smith   for (i=0; i<slice->Nghosts*bs; i++) {
7047c6ae99SBarry Smith     globals[slice->n*bs+i] = slice->ghosts[i/bs]*bs + i%bs;
7147c6ae99SBarry Smith   }
7247c6ae99SBarry Smith   ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_SELF,(slice->n+slice->Nghosts)*bs,globals,PETSC_OWN_POINTER,&lmap);CHKERRQ(ierr);
7347c6ae99SBarry Smith   ierr = ISLocalToGlobalMappingBlock(lmap,bs,&blmap);CHKERRQ(ierr);
74784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(*J,lmap,lmap);CHKERRQ(ierr);
75784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(*J,blmap,blmap);CHKERRQ(ierr);
76fcfd50ebSBarry Smith   ierr = ISLocalToGlobalMappingDestroy(&lmap);CHKERRQ(ierr);
77fcfd50ebSBarry Smith   ierr = ISLocalToGlobalMappingDestroy(&blmap);CHKERRQ(ierr);
7847c6ae99SBarry Smith   PetscFunctionReturn(0);
7947c6ae99SBarry Smith }
8047c6ae99SBarry Smith 
8147c6ae99SBarry Smith #undef __FUNCT__
820c010503SBarry Smith #define __FUNCT__ "DMSlicedSetGhosts"
8347c6ae99SBarry Smith /*@C
840c010503SBarry Smith     DMSlicedSetGhosts - Sets the global indices of other processes elements that will
8547c6ae99SBarry Smith       be ghosts on this process
8647c6ae99SBarry Smith 
8747c6ae99SBarry Smith     Not Collective
8847c6ae99SBarry Smith 
8947c6ae99SBarry Smith     Input Parameters:
900c010503SBarry Smith +    slice - the DM object
9147c6ae99SBarry Smith .    bs - block size
9247c6ae99SBarry Smith .    nlocal - number of local (owned, non-ghost) blocks
9347c6ae99SBarry Smith .    Nghosts - number of ghost blocks on this process
9447c6ae99SBarry Smith -    ghosts - global indices of each ghost block
9547c6ae99SBarry Smith 
9647c6ae99SBarry Smith     Level: advanced
9747c6ae99SBarry Smith 
98240edb01SJed Brown .seealso DMDestroy(), DMCreateGlobalVector()
9947c6ae99SBarry Smith 
10047c6ae99SBarry Smith @*/
1017087cfbeSBarry Smith PetscErrorCode  DMSlicedSetGhosts(DM dm,PetscInt bs,PetscInt nlocal,PetscInt Nghosts,const PetscInt ghosts[])
10247c6ae99SBarry Smith {
10347c6ae99SBarry Smith   PetscErrorCode ierr;
10447c6ae99SBarry Smith   DM_Sliced      *slice = (DM_Sliced*)dm->data;
10547c6ae99SBarry Smith 
10647c6ae99SBarry Smith   PetscFunctionBegin;
1070c010503SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
10847c6ae99SBarry Smith   ierr = PetscFree(slice->ghosts);CHKERRQ(ierr);
10947c6ae99SBarry Smith   ierr = PetscMalloc(Nghosts*sizeof(PetscInt),&slice->ghosts);CHKERRQ(ierr);
11047c6ae99SBarry Smith   ierr = PetscMemcpy(slice->ghosts,ghosts,Nghosts*sizeof(PetscInt));CHKERRQ(ierr);
11147c6ae99SBarry Smith   slice->bs      = bs;
11247c6ae99SBarry Smith   slice->n       = nlocal;
11347c6ae99SBarry Smith   slice->Nghosts = Nghosts;
11447c6ae99SBarry Smith   PetscFunctionReturn(0);
11547c6ae99SBarry Smith }
11647c6ae99SBarry Smith 
11747c6ae99SBarry Smith #undef __FUNCT__
1180c010503SBarry Smith #define __FUNCT__ "DMSlicedSetPreallocation"
11947c6ae99SBarry Smith /*@C
1200c010503SBarry Smith     DMSlicedSetPreallocation - sets the matrix memory preallocation for matrices computed by DMSliced
12147c6ae99SBarry Smith 
12247c6ae99SBarry Smith     Not Collective
12347c6ae99SBarry Smith 
12447c6ae99SBarry Smith     Input Parameters:
1250c010503SBarry Smith +    slice - the DM object
12647c6ae99SBarry Smith .    d_nz  - number of block nonzeros per block row in diagonal portion of local
12747c6ae99SBarry Smith            submatrix  (same for all local rows)
12847c6ae99SBarry Smith .    d_nnz - array containing the number of block nonzeros in the various block rows
12947c6ae99SBarry Smith            of the in diagonal portion of the local (possibly different for each block
13095742e49SBarry Smith            row) or PETSC_NULL.
13147c6ae99SBarry Smith .    o_nz  - number of block nonzeros per block row in the off-diagonal portion of local
13247c6ae99SBarry Smith            submatrix (same for all local rows).
13347c6ae99SBarry Smith -    o_nnz - array containing the number of nonzeros in the various block rows of the
13447c6ae99SBarry Smith            off-diagonal portion of the local submatrix (possibly different for
13547c6ae99SBarry Smith            each block row) or PETSC_NULL.
13647c6ae99SBarry Smith 
13747c6ae99SBarry Smith     Notes:
13847c6ae99SBarry Smith     See MatMPIBAIJSetPreallocation() for more details on preallocation.  If a scalar matrix (AIJ) is
1390c010503SBarry Smith     obtained with DMSlicedGetMatrix(), the correct preallocation will be set, respecting DMSlicedSetBlockFills().
14047c6ae99SBarry Smith 
14147c6ae99SBarry Smith     Level: advanced
14247c6ae99SBarry Smith 
143240edb01SJed Brown .seealso DMDestroy(), DMCreateGlobalVector(), MatMPIAIJSetPreallocation(),
1440c010503SBarry Smith          MatMPIBAIJSetPreallocation(), DMSlicedGetMatrix(), DMSlicedSetBlockFills()
14547c6ae99SBarry Smith 
14647c6ae99SBarry Smith @*/
1477087cfbeSBarry Smith PetscErrorCode  DMSlicedSetPreallocation(DM dm,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
14847c6ae99SBarry Smith {
14947c6ae99SBarry Smith   DM_Sliced *slice = (DM_Sliced*)dm->data;
15047c6ae99SBarry Smith 
15147c6ae99SBarry Smith   PetscFunctionBegin;
15247c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
15347c6ae99SBarry Smith   slice->d_nz  = d_nz;
15447c6ae99SBarry Smith   slice->d_nnz = (PetscInt*)d_nnz;
15547c6ae99SBarry Smith   slice->o_nz  = o_nz;
15647c6ae99SBarry Smith   slice->o_nnz = (PetscInt*)o_nnz;
15747c6ae99SBarry Smith   PetscFunctionReturn(0);
15847c6ae99SBarry Smith }
15947c6ae99SBarry Smith 
16047c6ae99SBarry Smith #undef __FUNCT__
1610c010503SBarry Smith #define __FUNCT__ "DMSlicedSetBlockFills_Private"
1620c010503SBarry Smith static PetscErrorCode DMSlicedSetBlockFills_Private(PetscInt bs,const PetscInt *fill,DMSlicedBlockFills **inf)
16347c6ae99SBarry Smith {
16447c6ae99SBarry Smith   PetscErrorCode     ierr;
16547c6ae99SBarry Smith   PetscInt           i,j,nz,*fi,*fj;
1660c010503SBarry Smith   DMSlicedBlockFills *f;
16747c6ae99SBarry Smith 
16847c6ae99SBarry Smith   PetscFunctionBegin;
16947c6ae99SBarry Smith   PetscValidPointer(inf,3);
17047c6ae99SBarry Smith   if (*inf) {ierr = PetscFree3((*inf)->i,(*inf)->j,*inf);CHKERRQ(ierr);}
17147c6ae99SBarry Smith   if (!fill) PetscFunctionReturn(0);
17247c6ae99SBarry Smith   for (i=0,nz=0; i<bs*bs; i++) if (fill[i]) nz++;
1730c010503SBarry Smith   ierr = PetscMalloc3(1,DMSlicedBlockFills,&f,bs+1,PetscInt,&fi,nz,PetscInt,&fj);CHKERRQ(ierr);
17447c6ae99SBarry Smith   f->bs = bs;
17547c6ae99SBarry Smith   f->nz = nz;
17647c6ae99SBarry Smith   f->i  = fi;
17747c6ae99SBarry Smith   f->j  = fj;
17847c6ae99SBarry Smith   for (i=0,nz=0; i<bs; i++) {
17947c6ae99SBarry Smith     fi[i] = nz;
18047c6ae99SBarry Smith     for (j=0; j<bs; j++) if (fill[i*bs+j]) fj[nz++] = j;
18147c6ae99SBarry Smith   }
18247c6ae99SBarry Smith   fi[i] = nz;
18347c6ae99SBarry Smith   *inf = f;
18447c6ae99SBarry Smith   PetscFunctionReturn(0);
18547c6ae99SBarry Smith }
18647c6ae99SBarry Smith 
18747c6ae99SBarry Smith #undef __FUNCT__
1880c010503SBarry Smith #define __FUNCT__ "DMSlicedSetBlockFills"
18947c6ae99SBarry Smith /*@C
1900c010503SBarry Smith     DMSlicedSetBlockFills - Sets the fill pattern in each block for a multi-component problem
1910c010503SBarry Smith     of the matrix returned by DMSlicedGetMatrix().
19247c6ae99SBarry Smith 
1930c010503SBarry Smith     Logically Collective on DM
19447c6ae99SBarry Smith 
19547c6ae99SBarry Smith     Input Parameter:
1960c010503SBarry Smith +   sliced - the DM object
19747c6ae99SBarry Smith .   dfill - the fill pattern in the diagonal block (may be PETSC_NULL, means use dense block)
19847c6ae99SBarry Smith -   ofill - the fill pattern in the off-diagonal blocks
19947c6ae99SBarry Smith 
20047c6ae99SBarry Smith     Notes:
20147c6ae99SBarry Smith     This only makes sense for multicomponent problems using scalar matrix formats (AIJ).
202aa219208SBarry Smith     See DMDASetBlockFills() for example usage.
20347c6ae99SBarry Smith 
20447c6ae99SBarry Smith     Level: advanced
20547c6ae99SBarry Smith 
206aa219208SBarry Smith .seealso DMSlicedGetMatrix(), DMDASetBlockFills()
20747c6ae99SBarry Smith @*/
2087087cfbeSBarry Smith PetscErrorCode  DMSlicedSetBlockFills(DM dm,const PetscInt *dfill,const PetscInt *ofill)
20947c6ae99SBarry Smith {
21047c6ae99SBarry Smith   DM_Sliced      *slice = (DM_Sliced*)dm->data;
21147c6ae99SBarry Smith   PetscErrorCode ierr;
21247c6ae99SBarry Smith 
21347c6ae99SBarry Smith   PetscFunctionBegin;
2140c010503SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2150c010503SBarry Smith   ierr = DMSlicedSetBlockFills_Private(slice->bs,dfill,&slice->dfill);CHKERRQ(ierr);
2160c010503SBarry Smith   ierr = DMSlicedSetBlockFills_Private(slice->bs,ofill,&slice->ofill);CHKERRQ(ierr);
21747c6ae99SBarry Smith   PetscFunctionReturn(0);
21847c6ae99SBarry Smith }
21947c6ae99SBarry Smith 
22047c6ae99SBarry Smith #undef __FUNCT__
2210c010503SBarry Smith #define __FUNCT__ "DMDestroy_Sliced"
2227087cfbeSBarry Smith PetscErrorCode  DMDestroy_Sliced(DM dm)
22347c6ae99SBarry Smith {
22447c6ae99SBarry Smith   PetscErrorCode ierr;
22547c6ae99SBarry Smith   DM_Sliced      *slice = (DM_Sliced*)dm->data;
22647c6ae99SBarry Smith 
22747c6ae99SBarry Smith   PetscFunctionBegin;
22847c6ae99SBarry Smith   ierr = PetscFree(slice->ghosts);CHKERRQ(ierr);
2293c4aec1bSBarry Smith   if (slice->dfill) {ierr = PetscFree3(slice->dfill,slice->dfill->i,slice->dfill->j);CHKERRQ(ierr);}
2303c4aec1bSBarry Smith   if (slice->ofill) {ierr = PetscFree3(slice->ofill,slice->ofill->i,slice->ofill->j);CHKERRQ(ierr);}
23147c6ae99SBarry Smith   PetscFunctionReturn(0);
23247c6ae99SBarry Smith }
23347c6ae99SBarry Smith 
23447c6ae99SBarry Smith #undef __FUNCT__
2350c010503SBarry Smith #define __FUNCT__ "DMCreateGlobalVector_Sliced"
2367087cfbeSBarry Smith PetscErrorCode  DMCreateGlobalVector_Sliced(DM dm,Vec *gvec)
23747c6ae99SBarry Smith {
23847c6ae99SBarry Smith   PetscErrorCode     ierr;
23947c6ae99SBarry Smith   DM_Sliced          *slice = (DM_Sliced*)dm->data;
24047c6ae99SBarry Smith 
24147c6ae99SBarry Smith   PetscFunctionBegin;
24247c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
24347c6ae99SBarry Smith   PetscValidPointer(gvec,2);
24447c6ae99SBarry Smith   *gvec = 0;
2459ff76d2cSJed Brown   ierr = VecCreateGhostBlock(((PetscObject)dm)->comm,slice->bs,slice->n*slice->bs,PETSC_DETERMINE,slice->Nghosts,slice->ghosts,gvec);CHKERRQ(ierr);
246b549a33dSJed Brown   ierr = PetscObjectCompose((PetscObject)*gvec,"DM",(PetscObject)dm);CHKERRQ(ierr);
24747c6ae99SBarry Smith   PetscFunctionReturn(0);
24847c6ae99SBarry Smith }
24947c6ae99SBarry Smith 
250a4121054SBarry Smith EXTERN_C_BEGIN
251a4121054SBarry Smith #undef __FUNCT__
252a4121054SBarry Smith #define __FUNCT__ "DMCreate_Sliced"
2537087cfbeSBarry Smith PetscErrorCode  DMCreate_Sliced(DM p)
254a4121054SBarry Smith {
255a4121054SBarry Smith   PetscErrorCode ierr;
256a4121054SBarry Smith   DM_Sliced      *slice;
257a4121054SBarry Smith 
258a4121054SBarry Smith   PetscFunctionBegin;
259a4121054SBarry Smith   ierr = PetscNewLog(p,DM_Sliced,&slice);CHKERRQ(ierr);
260a4121054SBarry Smith   p->data = slice;
261a4121054SBarry Smith 
262b549a33dSJed Brown   ierr = PetscObjectChangeTypeName((PetscObject)p,DMSLICED);CHKERRQ(ierr);
263a4121054SBarry Smith   p->ops->createglobalvector = DMCreateGlobalVector_Sliced;
26425296bd5SBarry Smith   p->ops->creatematrix          = DMCreateMatrix_Sliced;
265a4121054SBarry Smith   p->ops->destroy            = DMDestroy_Sliced;
266a4121054SBarry Smith   PetscFunctionReturn(0);
267a4121054SBarry Smith }
268a4121054SBarry Smith EXTERN_C_END
269a4121054SBarry Smith 
27047c6ae99SBarry Smith #undef __FUNCT__
2710c010503SBarry Smith #define __FUNCT__ "DMSlicedCreate"
27247c6ae99SBarry Smith /*@C
2730c010503SBarry Smith     DMSlicedCreate - Creates a DM object, used to manage data for a unstructured problem
27447c6ae99SBarry Smith 
2750c010503SBarry Smith     Collective on MPI_Comm
2760c010503SBarry Smith 
2770c010503SBarry Smith     Input Parameter:
2780c010503SBarry Smith .   comm - the processors that will share the global vector
2790c010503SBarry Smith 
2800c010503SBarry Smith     Output Parameters:
2810c010503SBarry Smith .   slice - the slice object
2820c010503SBarry Smith 
2830c010503SBarry Smith     Level: advanced
2840c010503SBarry Smith 
285240edb01SJed Brown .seealso DMDestroy(), DMCreateGlobalVector()
2860c010503SBarry Smith 
2870c010503SBarry Smith @*/
2887087cfbeSBarry Smith PetscErrorCode  DMSlicedCreate(MPI_Comm comm,DM *dm)
2890c010503SBarry Smith {
2900c010503SBarry Smith   PetscErrorCode ierr;
2910c010503SBarry Smith 
2920c010503SBarry Smith   PetscFunctionBegin;
2930c010503SBarry Smith   PetscValidPointer(dm,2);
294a4121054SBarry Smith   ierr = DMCreate(comm,dm);CHKERRQ(ierr);
295a4121054SBarry Smith   ierr = DMSetType(*dm,DMSLICED);CHKERRQ(ierr);
2960c010503SBarry Smith   PetscFunctionReturn(0);
2970c010503SBarry Smith }
2980c010503SBarry Smith 
299aa219208SBarry Smith /* Explanation of the missing functions for DMDA-style handling of the local vector:
30047c6ae99SBarry Smith 
3010c010503SBarry Smith    DMSlicedCreateLocalVector()
3020c010503SBarry Smith    DMSlicedGlobalToLocalBegin()
3030c010503SBarry Smith    DMSlicedGlobalToLocalEnd()
30447c6ae99SBarry Smith 
3050c010503SBarry Smith  There is no way to get the global form from a local form, so DMSlicedCreateLocalVector() is a memory leak without
3060c010503SBarry Smith  external accounting for the global vector.  Also, DMSliced intends the user to work with the VecGhost interface since the
307aa219208SBarry Smith  ghosts are already ordered after the owned entries.  Contrast this to a DMDA where the local vector has a special
30847c6ae99SBarry Smith  ordering described by the structured grid, hence it cannot share memory with the global form.  For this reason, users
3090c010503SBarry Smith  of DMSliced should work with the global vector and use
31047c6ae99SBarry Smith 
31147c6ae99SBarry Smith    VecGhostGetLocalForm(), VecGhostRestoreLocalForm()
31247c6ae99SBarry Smith    VecGhostUpdateBegin(), VecGhostUpdateEnd()
31347c6ae99SBarry Smith 
314aa219208SBarry Smith  rather than the missing DMDA-style functions.  This is conceptually simpler and offers better performance than is
315aa219208SBarry Smith  possible with the DMDA-style interface.
31647c6ae99SBarry Smith */
317