xref: /petsc/src/dm/impls/sliced/sliced.c (revision c688c0463f6c25494d0d243b2a43b2b3454e025c)
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"
1819fd82e9SBarry Smith PetscErrorCode  DMCreateMatrix_Sliced(DM dm, 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);
30f73d5cc4SBarry 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"
2223efe6655SBarry Smith static 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);}
231435a35e8SMatthew G Knepley   /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */
232435a35e8SMatthew G Knepley   ierr = PetscFree(slice);CHKERRQ(ierr);
23347c6ae99SBarry Smith   PetscFunctionReturn(0);
23447c6ae99SBarry Smith }
23547c6ae99SBarry Smith 
23647c6ae99SBarry Smith #undef __FUNCT__
2370c010503SBarry Smith #define __FUNCT__ "DMCreateGlobalVector_Sliced"
2383efe6655SBarry Smith static PetscErrorCode  DMCreateGlobalVector_Sliced(DM dm,Vec *gvec)
23947c6ae99SBarry Smith {
24047c6ae99SBarry Smith   PetscErrorCode     ierr;
24147c6ae99SBarry Smith   DM_Sliced          *slice = (DM_Sliced*)dm->data;
24247c6ae99SBarry Smith 
24347c6ae99SBarry Smith   PetscFunctionBegin;
24447c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
24547c6ae99SBarry Smith   PetscValidPointer(gvec,2);
24647c6ae99SBarry Smith   *gvec = 0;
2479ff76d2cSJed Brown   ierr = VecCreateGhostBlock(((PetscObject)dm)->comm,slice->bs,slice->n*slice->bs,PETSC_DETERMINE,slice->Nghosts,slice->ghosts,gvec);CHKERRQ(ierr);
248*c688c046SMatthew G Knepley   ierr = VecSetDM(*gvec,dm);CHKERRQ(ierr);
24947c6ae99SBarry Smith   PetscFunctionReturn(0);
25047c6ae99SBarry Smith }
25147c6ae99SBarry Smith 
2523efe6655SBarry Smith #undef __FUNCT__
2533efe6655SBarry Smith #define __FUNCT__ "DMGlobalToLocalBegin_Sliced"
2543efe6655SBarry Smith static PetscErrorCode  DMGlobalToLocalBegin_Sliced(DM da,Vec g,InsertMode mode,Vec l)
2553efe6655SBarry Smith {
2563efe6655SBarry Smith   PetscErrorCode ierr;
2573efe6655SBarry Smith   PetscBool      flg;
2583efe6655SBarry Smith 
2593efe6655SBarry Smith   PetscFunctionBegin;
260a3314f2cSMatthew G Knepley   ierr = VecGhostIsLocalForm(g,l,&flg);CHKERRQ(ierr);
2613efe6655SBarry Smith   if (!flg) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Local vector is not local form of global vector");
2623efe6655SBarry Smith   ierr = VecGhostUpdateEnd(g,mode,SCATTER_FORWARD);CHKERRQ(ierr);
2633efe6655SBarry Smith   ierr = VecGhostUpdateBegin(g,mode,SCATTER_FORWARD);CHKERRQ(ierr);
2643efe6655SBarry Smith   PetscFunctionReturn(0);
2653efe6655SBarry Smith }
2663efe6655SBarry Smith 
2673efe6655SBarry Smith #undef __FUNCT__
2683efe6655SBarry Smith #define __FUNCT__ "DMGlobalToLocalEnd_Sliced"
2693efe6655SBarry Smith static PetscErrorCode  DMGlobalToLocalEnd_Sliced(DM da,Vec g,InsertMode mode,Vec l)
2703efe6655SBarry Smith {
2713efe6655SBarry Smith   PetscErrorCode ierr;
2723efe6655SBarry Smith   PetscBool      flg;
2733efe6655SBarry Smith 
2743efe6655SBarry Smith   PetscFunctionBegin;
2753efe6655SBarry Smith   ierr = VecGhostIsLocalForm(g,l,&flg);CHKERRQ(ierr);
2763efe6655SBarry Smith   if (!flg) SETERRQ(((PetscObject)da)->comm,PETSC_ERR_ARG_WRONG,"Local vector is not local form of global vector");
2773efe6655SBarry Smith   ierr = VecGhostUpdateEnd(g,mode,SCATTER_FORWARD);CHKERRQ(ierr);
2783efe6655SBarry Smith   PetscFunctionReturn(0);
2793efe6655SBarry Smith }
2803efe6655SBarry Smith 
2813e0de288SBarry Smith /*MC
2823e0de288SBarry Smith    DMSLICED = "sliced" - A DM object that is used to manage data for a general graph. Uses VecCreateGhost() ghosted vectors for storing the fields
2833e0de288SBarry Smith 
2843e0de288SBarry Smith    See DMCreateSliced() for details.
2853e0de288SBarry Smith 
2863e0de288SBarry Smith   Level: intermediate
2873e0de288SBarry Smith 
2883e0de288SBarry Smith .seealso: DMType, DMCOMPOSITE, DMCreateSliced(), DMCreate()
2893e0de288SBarry Smith M*/
2903e0de288SBarry Smith 
291a4121054SBarry Smith EXTERN_C_BEGIN
292a4121054SBarry Smith #undef __FUNCT__
293a4121054SBarry Smith #define __FUNCT__ "DMCreate_Sliced"
2947087cfbeSBarry Smith PetscErrorCode  DMCreate_Sliced(DM p)
295a4121054SBarry Smith {
296a4121054SBarry Smith   PetscErrorCode ierr;
297a4121054SBarry Smith   DM_Sliced      *slice;
298a4121054SBarry Smith 
299a4121054SBarry Smith   PetscFunctionBegin;
300a4121054SBarry Smith   ierr = PetscNewLog(p,DM_Sliced,&slice);CHKERRQ(ierr);
301a4121054SBarry Smith   p->data = slice;
302a4121054SBarry Smith 
303b549a33dSJed Brown   ierr = PetscObjectChangeTypeName((PetscObject)p,DMSLICED);CHKERRQ(ierr);
304a4121054SBarry Smith   p->ops->createglobalvector = DMCreateGlobalVector_Sliced;
30525296bd5SBarry Smith   p->ops->creatematrix       = DMCreateMatrix_Sliced;
3063efe6655SBarry Smith   p->ops->globaltolocalbegin = DMGlobalToLocalBegin_Sliced;
3073efe6655SBarry Smith   p->ops->globaltolocalend   = DMGlobalToLocalEnd_Sliced;
308a4121054SBarry Smith   p->ops->destroy            = DMDestroy_Sliced;
309a4121054SBarry Smith   PetscFunctionReturn(0);
310a4121054SBarry Smith }
311a4121054SBarry Smith EXTERN_C_END
312a4121054SBarry Smith 
31347c6ae99SBarry Smith #undef __FUNCT__
3140c010503SBarry Smith #define __FUNCT__ "DMSlicedCreate"
31547c6ae99SBarry Smith /*@C
3160c010503SBarry Smith     DMSlicedCreate - Creates a DM object, used to manage data for a unstructured problem
31747c6ae99SBarry Smith 
3180c010503SBarry Smith     Collective on MPI_Comm
3190c010503SBarry Smith 
3200c010503SBarry Smith     Input Parameter:
3213efe6655SBarry Smith +   comm - the processors that will share the global vector
3223efe6655SBarry Smith .   bs - the block size
3233efe6655SBarry Smith .   nlocal - number of vector entries on this process
3243efe6655SBarry Smith .   Nghosts - number of ghost points needed on this process
3253efe6655SBarry Smith .   ghosts - global indices of all ghost points for this process
3263efe6655SBarry Smith .   d_nnz - matrix preallocation information representing coupling within this process
3273efe6655SBarry Smith -   o_nnz - matrix preallocation information representing coupling between this process and other processes
3280c010503SBarry Smith 
3290c010503SBarry Smith     Output Parameters:
3300c010503SBarry Smith .   slice - the slice object
3310c010503SBarry Smith 
3323efe6655SBarry Smith     Notes:
3333efe6655SBarry Smith         This DM does not support DMCreateLocalVector(), DMGlobalToLocalBegin(), and DMGlobalToLocalEnd() instead one directly uses
3343efe6655SBarry Smith         VecGhostGetLocalForm() and VecGhostRestoreLocalForm() to access the local representation and VecGhostUpdateBegin() and VecGhostUpdateEnd() to update
3353efe6655SBarry Smith         the ghost points.
3363efe6655SBarry Smith 
3373efe6655SBarry Smith         One can use DMGlobalToLocalBegin(), and DMGlobalToLocalEnd() instead of VecGhostUpdateBegin() and VecGhostUpdateEnd().
3383efe6655SBarry Smith 
3390c010503SBarry Smith     Level: advanced
3400c010503SBarry Smith 
3413e0de288SBarry Smith .seealso DMDestroy(), DMCreateGlobalVector(), DMSetType(), DMSLICED, DMSlicedSetGhosts(), DMSlicedSetPreallocation(), VecGhostUpdateBegin(), VecGhostUpdateEnd(),
3423e0de288SBarry Smith          VecGhostGetLocalForm(), VecGhostRestoreLocalForm()
3430c010503SBarry Smith 
3440c010503SBarry Smith @*/
3453efe6655SBarry 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)
3460c010503SBarry Smith {
3470c010503SBarry Smith   PetscErrorCode ierr;
3480c010503SBarry Smith 
3490c010503SBarry Smith   PetscFunctionBegin;
3500c010503SBarry Smith   PetscValidPointer(dm,2);
351a4121054SBarry Smith   ierr = DMCreate(comm,dm);CHKERRQ(ierr);
352a4121054SBarry Smith   ierr = DMSetType(*dm,DMSLICED);CHKERRQ(ierr);
3533efe6655SBarry Smith   ierr = DMSlicedSetGhosts(*dm,bs,nlocal,Nghosts,ghosts);CHKERRQ(ierr);
3543efe6655SBarry Smith   if (d_nnz) {
3553efe6655SBarry Smith     ierr = DMSlicedSetPreallocation(*dm,0, d_nnz,0,o_nnz);CHKERRQ(ierr);
3563efe6655SBarry Smith   }
3570c010503SBarry Smith   PetscFunctionReturn(0);
3580c010503SBarry Smith }
3590c010503SBarry Smith 
360