xref: /petsc/src/dm/impls/sliced/sliced.c (revision b549a33df52e18f44cd122c2139f7c5486754f82)
1c6db04a5SJed Brown #include <petscdmsliced.h>      /*I      "petscdmsliced.h"     I*/
2c6db04a5SJed Brown #include <petscmat.h>           /*I      "petscmat.h"    I*/
3c6db04a5SJed Brown #include <private/dmimpl.h>     /*I      "petscmat.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   Vec                globalvector;
1247c6ae99SBarry Smith   PetscInt           bs,n,N,Nghosts,*ghosts;
1347c6ae99SBarry Smith   PetscInt           d_nz,o_nz,*d_nnz,*o_nnz;
140c010503SBarry Smith   DMSlicedBlockFills *dfill,*ofill;
1547c6ae99SBarry Smith } DM_Sliced;
1647c6ae99SBarry Smith 
1747c6ae99SBarry Smith #undef __FUNCT__
180c010503SBarry Smith #define __FUNCT__ "DMGetMatrix_Sliced"
197087cfbeSBarry Smith PetscErrorCode  DMGetMatrix_Sliced(DM dm, const MatType mtype,Mat *J)
2047c6ae99SBarry Smith {
2147c6ae99SBarry Smith   PetscErrorCode         ierr;
2247c6ae99SBarry Smith   PetscInt               *globals,*sd_nnz,*so_nnz,rstart,bs,i;
2347c6ae99SBarry Smith   ISLocalToGlobalMapping lmap,blmap;
2447c6ae99SBarry Smith   void                   (*aij)(void) = PETSC_NULL;
2547c6ae99SBarry Smith   DM_Sliced              *slice = (DM_Sliced*)dm->data;
2647c6ae99SBarry Smith 
2747c6ae99SBarry Smith   PetscFunctionBegin;
2847c6ae99SBarry Smith   bs = slice->bs;
2947c6ae99SBarry Smith   ierr = MatCreate(((PetscObject)dm)->comm,J);CHKERRQ(ierr);
3047c6ae99SBarry Smith   ierr = MatSetSizes(*J,slice->n*bs,slice->n*bs,PETSC_DETERMINE,PETSC_DETERMINE);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   ierr = MatSetBlockSize(*J,bs);CHKERRQ(ierr);
6447c6ae99SBarry Smith 
6547c6ae99SBarry Smith   /* Set up the local to global map.  For the scalar map, we have to translate to entry-wise indexing instead of block-wise. */
6647c6ae99SBarry Smith   ierr = PetscMalloc((slice->n+slice->Nghosts)*bs*sizeof(PetscInt),&globals);CHKERRQ(ierr);
6747c6ae99SBarry Smith   ierr = MatGetOwnershipRange(*J,&rstart,PETSC_NULL);CHKERRQ(ierr);
6847c6ae99SBarry Smith   for (i=0; i<slice->n*bs; i++) {
6947c6ae99SBarry Smith     globals[i] = rstart + i;
7047c6ae99SBarry Smith   }
7147c6ae99SBarry Smith   for (i=0; i<slice->Nghosts*bs; i++) {
7247c6ae99SBarry Smith     globals[slice->n*bs+i] = slice->ghosts[i/bs]*bs + i%bs;
7347c6ae99SBarry Smith   }
7447c6ae99SBarry Smith   ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_SELF,(slice->n+slice->Nghosts)*bs,globals,PETSC_OWN_POINTER,&lmap);CHKERRQ(ierr);
7547c6ae99SBarry Smith   ierr = ISLocalToGlobalMappingBlock(lmap,bs,&blmap);CHKERRQ(ierr);
76784ac674SJed Brown   ierr = MatSetLocalToGlobalMapping(*J,lmap,lmap);CHKERRQ(ierr);
77784ac674SJed Brown   ierr = MatSetLocalToGlobalMappingBlock(*J,blmap,blmap);CHKERRQ(ierr);
78fcfd50ebSBarry Smith   ierr = ISLocalToGlobalMappingDestroy(&lmap);CHKERRQ(ierr);
79fcfd50ebSBarry Smith   ierr = ISLocalToGlobalMappingDestroy(&blmap);CHKERRQ(ierr);
8047c6ae99SBarry Smith   PetscFunctionReturn(0);
8147c6ae99SBarry Smith }
8247c6ae99SBarry Smith 
8347c6ae99SBarry Smith #undef __FUNCT__
840c010503SBarry Smith #define __FUNCT__ "DMSlicedSetGhosts"
8547c6ae99SBarry Smith /*@C
860c010503SBarry Smith     DMSlicedSetGhosts - Sets the global indices of other processes elements that will
8747c6ae99SBarry Smith       be ghosts on this process
8847c6ae99SBarry Smith 
8947c6ae99SBarry Smith     Not Collective
9047c6ae99SBarry Smith 
9147c6ae99SBarry Smith     Input Parameters:
920c010503SBarry Smith +    slice - the DM object
9347c6ae99SBarry Smith .    bs - block size
9447c6ae99SBarry Smith .    nlocal - number of local (owned, non-ghost) blocks
9547c6ae99SBarry Smith .    Nghosts - number of ghost blocks on this process
9647c6ae99SBarry Smith -    ghosts - global indices of each ghost block
9747c6ae99SBarry Smith 
9847c6ae99SBarry Smith     Level: advanced
9947c6ae99SBarry Smith 
1000c010503SBarry Smith .seealso DMDestroy(), DMCreateGlobalVector(), DMSlicedGetGlobalIndices()
10147c6ae99SBarry Smith 
10247c6ae99SBarry Smith @*/
1037087cfbeSBarry Smith PetscErrorCode  DMSlicedSetGhosts(DM dm,PetscInt bs,PetscInt nlocal,PetscInt Nghosts,const PetscInt ghosts[])
10447c6ae99SBarry Smith {
10547c6ae99SBarry Smith   PetscErrorCode ierr;
10647c6ae99SBarry Smith   DM_Sliced      *slice = (DM_Sliced*)dm->data;
10747c6ae99SBarry Smith 
10847c6ae99SBarry Smith   PetscFunctionBegin;
1090c010503SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
11047c6ae99SBarry Smith   ierr = PetscFree(slice->ghosts);CHKERRQ(ierr);
11147c6ae99SBarry Smith   ierr = PetscMalloc(Nghosts*sizeof(PetscInt),&slice->ghosts);CHKERRQ(ierr);
11247c6ae99SBarry Smith   ierr = PetscMemcpy(slice->ghosts,ghosts,Nghosts*sizeof(PetscInt));CHKERRQ(ierr);
11347c6ae99SBarry Smith   slice->bs      = bs;
11447c6ae99SBarry Smith   slice->n       = nlocal;
11547c6ae99SBarry Smith   slice->Nghosts = Nghosts;
11647c6ae99SBarry Smith   PetscFunctionReturn(0);
11747c6ae99SBarry Smith }
11847c6ae99SBarry Smith 
11947c6ae99SBarry Smith #undef __FUNCT__
1200c010503SBarry Smith #define __FUNCT__ "DMSlicedSetPreallocation"
12147c6ae99SBarry Smith /*@C
1220c010503SBarry Smith     DMSlicedSetPreallocation - sets the matrix memory preallocation for matrices computed by DMSliced
12347c6ae99SBarry Smith 
12447c6ae99SBarry Smith     Not Collective
12547c6ae99SBarry Smith 
12647c6ae99SBarry Smith     Input Parameters:
1270c010503SBarry Smith +    slice - the DM object
12847c6ae99SBarry Smith .    d_nz  - number of block nonzeros per block row in diagonal portion of local
12947c6ae99SBarry Smith            submatrix  (same for all local rows)
13047c6ae99SBarry Smith .    d_nnz - array containing the number of block nonzeros in the various block rows
13147c6ae99SBarry Smith            of the in diagonal portion of the local (possibly different for each block
13247c6ae99SBarry Smith            row) or PETSC_NULL.  You must leave room for the diagonal entry even if it is zero.
13347c6ae99SBarry Smith .    o_nz  - number of block nonzeros per block row in the off-diagonal portion of local
13447c6ae99SBarry Smith            submatrix (same for all local rows).
13547c6ae99SBarry Smith -    o_nnz - array containing the number of nonzeros in the various block rows of the
13647c6ae99SBarry Smith            off-diagonal portion of the local submatrix (possibly different for
13747c6ae99SBarry Smith            each block row) or PETSC_NULL.
13847c6ae99SBarry Smith 
13947c6ae99SBarry Smith     Notes:
14047c6ae99SBarry Smith     See MatMPIBAIJSetPreallocation() for more details on preallocation.  If a scalar matrix (AIJ) is
1410c010503SBarry Smith     obtained with DMSlicedGetMatrix(), the correct preallocation will be set, respecting DMSlicedSetBlockFills().
14247c6ae99SBarry Smith 
14347c6ae99SBarry Smith     Level: advanced
14447c6ae99SBarry Smith 
1450c010503SBarry Smith .seealso DMDestroy(), DMCreateGlobalVector(), DMSlicedGetGlobalIndices(), MatMPIAIJSetPreallocation(),
1460c010503SBarry Smith          MatMPIBAIJSetPreallocation(), DMSlicedGetMatrix(), DMSlicedSetBlockFills()
14747c6ae99SBarry Smith 
14847c6ae99SBarry Smith @*/
1497087cfbeSBarry Smith PetscErrorCode  DMSlicedSetPreallocation(DM dm,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[])
15047c6ae99SBarry Smith {
15147c6ae99SBarry Smith   DM_Sliced *slice = (DM_Sliced*)dm->data;
15247c6ae99SBarry Smith 
15347c6ae99SBarry Smith   PetscFunctionBegin;
15447c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
15547c6ae99SBarry Smith   slice->d_nz  = d_nz;
15647c6ae99SBarry Smith   slice->d_nnz = (PetscInt*)d_nnz;
15747c6ae99SBarry Smith   slice->o_nz  = o_nz;
15847c6ae99SBarry Smith   slice->o_nnz = (PetscInt*)o_nnz;
15947c6ae99SBarry Smith   PetscFunctionReturn(0);
16047c6ae99SBarry Smith }
16147c6ae99SBarry Smith 
16247c6ae99SBarry Smith #undef __FUNCT__
1630c010503SBarry Smith #define __FUNCT__ "DMSlicedSetBlockFills_Private"
1640c010503SBarry Smith static PetscErrorCode DMSlicedSetBlockFills_Private(PetscInt bs,const PetscInt *fill,DMSlicedBlockFills **inf)
16547c6ae99SBarry Smith {
16647c6ae99SBarry Smith   PetscErrorCode     ierr;
16747c6ae99SBarry Smith   PetscInt           i,j,nz,*fi,*fj;
1680c010503SBarry Smith   DMSlicedBlockFills *f;
16947c6ae99SBarry Smith 
17047c6ae99SBarry Smith   PetscFunctionBegin;
17147c6ae99SBarry Smith   PetscValidPointer(inf,3);
17247c6ae99SBarry Smith   if (*inf) {ierr = PetscFree3((*inf)->i,(*inf)->j,*inf);CHKERRQ(ierr);}
17347c6ae99SBarry Smith   if (!fill) PetscFunctionReturn(0);
17447c6ae99SBarry Smith   for (i=0,nz=0; i<bs*bs; i++) if (fill[i]) nz++;
1750c010503SBarry Smith   ierr = PetscMalloc3(1,DMSlicedBlockFills,&f,bs+1,PetscInt,&fi,nz,PetscInt,&fj);CHKERRQ(ierr);
17647c6ae99SBarry Smith   f->bs = bs;
17747c6ae99SBarry Smith   f->nz = nz;
17847c6ae99SBarry Smith   f->i  = fi;
17947c6ae99SBarry Smith   f->j  = fj;
18047c6ae99SBarry Smith   for (i=0,nz=0; i<bs; i++) {
18147c6ae99SBarry Smith     fi[i] = nz;
18247c6ae99SBarry Smith     for (j=0; j<bs; j++) if (fill[i*bs+j]) fj[nz++] = j;
18347c6ae99SBarry Smith   }
18447c6ae99SBarry Smith   fi[i] = nz;
18547c6ae99SBarry Smith   *inf = f;
18647c6ae99SBarry Smith   PetscFunctionReturn(0);
18747c6ae99SBarry Smith }
18847c6ae99SBarry Smith 
18947c6ae99SBarry Smith #undef __FUNCT__
1900c010503SBarry Smith #define __FUNCT__ "DMSlicedSetBlockFills"
19147c6ae99SBarry Smith /*@C
1920c010503SBarry Smith     DMSlicedSetBlockFills - Sets the fill pattern in each block for a multi-component problem
1930c010503SBarry Smith     of the matrix returned by DMSlicedGetMatrix().
19447c6ae99SBarry Smith 
1950c010503SBarry Smith     Logically Collective on DM
19647c6ae99SBarry Smith 
19747c6ae99SBarry Smith     Input Parameter:
1980c010503SBarry Smith +   sliced - the DM object
19947c6ae99SBarry Smith .   dfill - the fill pattern in the diagonal block (may be PETSC_NULL, means use dense block)
20047c6ae99SBarry Smith -   ofill - the fill pattern in the off-diagonal blocks
20147c6ae99SBarry Smith 
20247c6ae99SBarry Smith     Notes:
20347c6ae99SBarry Smith     This only makes sense for multicomponent problems using scalar matrix formats (AIJ).
204aa219208SBarry Smith     See DMDASetBlockFills() for example usage.
20547c6ae99SBarry Smith 
20647c6ae99SBarry Smith     Level: advanced
20747c6ae99SBarry Smith 
208aa219208SBarry Smith .seealso DMSlicedGetMatrix(), DMDASetBlockFills()
20947c6ae99SBarry Smith @*/
2107087cfbeSBarry Smith PetscErrorCode  DMSlicedSetBlockFills(DM dm,const PetscInt *dfill,const PetscInt *ofill)
21147c6ae99SBarry Smith {
21247c6ae99SBarry Smith   DM_Sliced      *slice = (DM_Sliced*)dm->data;
21347c6ae99SBarry Smith   PetscErrorCode ierr;
21447c6ae99SBarry Smith 
21547c6ae99SBarry Smith   PetscFunctionBegin;
2160c010503SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
2170c010503SBarry Smith   ierr = DMSlicedSetBlockFills_Private(slice->bs,dfill,&slice->dfill);CHKERRQ(ierr);
2180c010503SBarry Smith   ierr = DMSlicedSetBlockFills_Private(slice->bs,ofill,&slice->ofill);CHKERRQ(ierr);
21947c6ae99SBarry Smith   PetscFunctionReturn(0);
22047c6ae99SBarry Smith }
22147c6ae99SBarry Smith 
22247c6ae99SBarry Smith #undef __FUNCT__
2230c010503SBarry Smith #define __FUNCT__ "DMDestroy_Sliced"
2247087cfbeSBarry Smith PetscErrorCode  DMDestroy_Sliced(DM dm)
22547c6ae99SBarry Smith {
22647c6ae99SBarry Smith   PetscErrorCode ierr;
22747c6ae99SBarry Smith   DM_Sliced      *slice = (DM_Sliced*)dm->data;
22847c6ae99SBarry Smith 
22947c6ae99SBarry Smith   PetscFunctionBegin;
2306bf464f9SBarry Smith   ierr = VecDestroy(&slice->globalvector);CHKERRQ(ierr);
23147c6ae99SBarry Smith   ierr = PetscFree(slice->ghosts);CHKERRQ(ierr);
2323c4aec1bSBarry Smith   if (slice->dfill) {ierr = PetscFree3(slice->dfill,slice->dfill->i,slice->dfill->j);CHKERRQ(ierr);}
2333c4aec1bSBarry Smith   if (slice->ofill) {ierr = PetscFree3(slice->ofill,slice->ofill->i,slice->ofill->j);CHKERRQ(ierr);}
23447c6ae99SBarry Smith   PetscFunctionReturn(0);
23547c6ae99SBarry Smith }
23647c6ae99SBarry Smith 
23747c6ae99SBarry Smith #undef __FUNCT__
2380c010503SBarry Smith #define __FUNCT__ "DMCreateGlobalVector_Sliced"
2397087cfbeSBarry Smith PetscErrorCode  DMCreateGlobalVector_Sliced(DM dm,Vec *gvec)
24047c6ae99SBarry Smith {
24147c6ae99SBarry Smith   PetscErrorCode     ierr;
24247c6ae99SBarry Smith   PetscInt           bs,cnt;
24347c6ae99SBarry Smith   DM_Sliced          *slice = (DM_Sliced*)dm->data;
24447c6ae99SBarry Smith 
24547c6ae99SBarry Smith   PetscFunctionBegin;
24647c6ae99SBarry Smith   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
24747c6ae99SBarry Smith   PetscValidPointer(gvec,2);
24847c6ae99SBarry Smith   *gvec = 0;
24947c6ae99SBarry Smith   if (slice->globalvector) {
25047c6ae99SBarry Smith     ierr = PetscObjectGetReference((PetscObject)slice->globalvector,&cnt);CHKERRQ(ierr);
25147c6ae99SBarry Smith     if (cnt == 1) {             /* Nobody else has a reference so we can just reference it and give it away */
25247c6ae99SBarry Smith       *gvec = slice->globalvector;
25347c6ae99SBarry Smith       ierr = PetscObjectReference((PetscObject)*gvec);CHKERRQ(ierr);
25447c6ae99SBarry Smith       ierr = VecZeroEntries(*gvec);CHKERRQ(ierr);
25547c6ae99SBarry Smith     } else {                    /* Someone else has a reference so we duplicate the global vector */
25647c6ae99SBarry Smith       ierr = VecDuplicate(slice->globalvector,gvec);CHKERRQ(ierr);
25747c6ae99SBarry Smith     }
25847c6ae99SBarry Smith   } else {
25947c6ae99SBarry Smith     bs = slice->bs;
26047c6ae99SBarry Smith     ierr = VecCreateGhostBlock(((PetscObject)dm)->comm,bs,slice->n*bs,PETSC_DETERMINE,slice->Nghosts,slice->ghosts,&slice->globalvector);CHKERRQ(ierr);
26147c6ae99SBarry Smith     *gvec = slice->globalvector;
262*b549a33dSJed Brown     ierr = PetscObjectCompose((PetscObject)*gvec,"DM",(PetscObject)dm);CHKERRQ(ierr);
26347c6ae99SBarry Smith     ierr = PetscObjectReference((PetscObject)*gvec);CHKERRQ(ierr);
26447c6ae99SBarry Smith   }
26547c6ae99SBarry Smith   PetscFunctionReturn(0);
26647c6ae99SBarry Smith }
26747c6ae99SBarry Smith 
268a4121054SBarry Smith EXTERN_C_BEGIN
269a4121054SBarry Smith #undef __FUNCT__
270a4121054SBarry Smith #define __FUNCT__ "DMCreate_Sliced"
2717087cfbeSBarry Smith PetscErrorCode  DMCreate_Sliced(DM p)
272a4121054SBarry Smith {
273a4121054SBarry Smith   PetscErrorCode ierr;
274a4121054SBarry Smith   DM_Sliced      *slice;
275a4121054SBarry Smith 
276a4121054SBarry Smith   PetscFunctionBegin;
277a4121054SBarry Smith   ierr = PetscNewLog(p,DM_Sliced,&slice);CHKERRQ(ierr);
278a4121054SBarry Smith   p->data = slice;
279a4121054SBarry Smith 
280*b549a33dSJed Brown   ierr = PetscObjectChangeTypeName((PetscObject)p,DMSLICED);CHKERRQ(ierr);
281a4121054SBarry Smith   p->ops->createglobalvector = DMCreateGlobalVector_Sliced;
282a4121054SBarry Smith   p->ops->getmatrix          = DMGetMatrix_Sliced;
283a4121054SBarry Smith   p->ops->destroy            = DMDestroy_Sliced;
284a4121054SBarry Smith   PetscFunctionReturn(0);
285a4121054SBarry Smith }
286a4121054SBarry Smith EXTERN_C_END
287a4121054SBarry Smith 
28847c6ae99SBarry Smith #undef __FUNCT__
2890c010503SBarry Smith #define __FUNCT__ "DMSlicedCreate"
29047c6ae99SBarry Smith /*@C
2910c010503SBarry Smith     DMSlicedCreate - Creates a DM object, used to manage data for a unstructured problem
29247c6ae99SBarry Smith 
2930c010503SBarry Smith     Collective on MPI_Comm
2940c010503SBarry Smith 
2950c010503SBarry Smith     Input Parameter:
2960c010503SBarry Smith .   comm - the processors that will share the global vector
2970c010503SBarry Smith 
2980c010503SBarry Smith     Output Parameters:
2990c010503SBarry Smith .   slice - the slice object
3000c010503SBarry Smith 
3010c010503SBarry Smith     Level: advanced
3020c010503SBarry Smith 
3030c010503SBarry Smith .seealso DMDestroy(), DMCreateGlobalVector(), DMSlicedGetGlobalIndices()
3040c010503SBarry Smith 
3050c010503SBarry Smith @*/
3067087cfbeSBarry Smith PetscErrorCode  DMSlicedCreate(MPI_Comm comm,DM *dm)
3070c010503SBarry Smith {
3080c010503SBarry Smith   PetscErrorCode ierr;
3090c010503SBarry Smith 
3100c010503SBarry Smith   PetscFunctionBegin;
3110c010503SBarry Smith   PetscValidPointer(dm,2);
312a4121054SBarry Smith   ierr = DMCreate(comm,dm);CHKERRQ(ierr);
313a4121054SBarry Smith   ierr = DMSetType(*dm,DMSLICED);CHKERRQ(ierr);
3140c010503SBarry Smith   PetscFunctionReturn(0);
3150c010503SBarry Smith }
3160c010503SBarry Smith 
3170c010503SBarry Smith 
3180c010503SBarry Smith #undef __FUNCT__
3190c010503SBarry Smith #define __FUNCT__ "DMSlicedGetGlobalIndices"
3200c010503SBarry Smith /*@C
3210c010503SBarry Smith     DMSlicedGetGlobalIndices - Gets the global indices for all the local entries
3220c010503SBarry Smith 
3230c010503SBarry Smith     Collective on DM
32447c6ae99SBarry Smith 
32547c6ae99SBarry Smith     Input Parameter:
32647c6ae99SBarry Smith .    slice - the slice object
32747c6ae99SBarry Smith 
32847c6ae99SBarry Smith     Output Parameters:
32947c6ae99SBarry Smith .    idx - the individual indices for each packed vector/array
33047c6ae99SBarry Smith 
33147c6ae99SBarry Smith     Level: advanced
33247c6ae99SBarry Smith 
33347c6ae99SBarry Smith     Notes:
33447c6ae99SBarry Smith        The idx parameters should be freed by the calling routine with PetscFree()
33547c6ae99SBarry Smith 
336*b549a33dSJed Brown .seealso DMSlicedDestroy(), DMCreateGlobalVector(), DMSlicedCreate()
33747c6ae99SBarry Smith 
33847c6ae99SBarry Smith @*/
3397087cfbeSBarry Smith PetscErrorCode  DMSlicedGetGlobalIndices(DM dm,PetscInt *idx[])
34047c6ae99SBarry Smith {
34147c6ae99SBarry Smith   PetscFunctionReturn(0);
34247c6ae99SBarry Smith }
34347c6ae99SBarry Smith 
34447c6ae99SBarry Smith 
345aa219208SBarry Smith /* Explanation of the missing functions for DMDA-style handling of the local vector:
34647c6ae99SBarry Smith 
3470c010503SBarry Smith    DMSlicedCreateLocalVector()
3480c010503SBarry Smith    DMSlicedGlobalToLocalBegin()
3490c010503SBarry Smith    DMSlicedGlobalToLocalEnd()
35047c6ae99SBarry Smith 
3510c010503SBarry Smith  There is no way to get the global form from a local form, so DMSlicedCreateLocalVector() is a memory leak without
3520c010503SBarry Smith  external accounting for the global vector.  Also, DMSliced intends the user to work with the VecGhost interface since the
353aa219208SBarry Smith  ghosts are already ordered after the owned entries.  Contrast this to a DMDA where the local vector has a special
35447c6ae99SBarry Smith  ordering described by the structured grid, hence it cannot share memory with the global form.  For this reason, users
3550c010503SBarry Smith  of DMSliced should work with the global vector and use
35647c6ae99SBarry Smith 
35747c6ae99SBarry Smith    VecGhostGetLocalForm(), VecGhostRestoreLocalForm()
35847c6ae99SBarry Smith    VecGhostUpdateBegin(), VecGhostUpdateEnd()
35947c6ae99SBarry Smith 
360aa219208SBarry Smith  rather than the missing DMDA-style functions.  This is conceptually simpler and offers better performance than is
361aa219208SBarry Smith  possible with the DMDA-style interface.
36247c6ae99SBarry Smith */
363