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 PetscErrorCode ierr; 1947c6ae99SBarry Smith PetscInt *globals,*sd_nnz,*so_nnz,rstart,bs,i; 2045b6f7e9SBarry Smith ISLocalToGlobalMapping lmap; 21f77a5242SKarl Rupp void (*aij)(void) = NULL; 2247c6ae99SBarry Smith DM_Sliced *slice = (DM_Sliced*)dm->data; 2347c6ae99SBarry Smith 2447c6ae99SBarry Smith PetscFunctionBegin; 2547c6ae99SBarry Smith bs = slice->bs; 26ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)dm),J);CHKERRQ(ierr); 2747c6ae99SBarry Smith ierr = MatSetSizes(*J,slice->n*bs,slice->n*bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 28f73d5cc4SBarry Smith ierr = MatSetBlockSize(*J,bs);CHKERRQ(ierr); 29b412c318SBarry Smith ierr = MatSetType(*J,dm->mattype);CHKERRQ(ierr); 3047c6ae99SBarry Smith ierr = MatSeqBAIJSetPreallocation(*J,bs,slice->d_nz,slice->d_nnz);CHKERRQ(ierr); 3147c6ae99SBarry Smith ierr = MatMPIBAIJSetPreallocation(*J,bs,slice->d_nz,slice->d_nnz,slice->o_nz,slice->o_nnz);CHKERRQ(ierr); 3247c6ae99SBarry Smith /* In general, we have to do extra work to preallocate for scalar (AIJ) matrices so we check whether it will do any 3347c6ae99SBarry Smith * good before going on with it. */ 3447c6ae99SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)*J,"MatMPIAIJSetPreallocation_C",&aij);CHKERRQ(ierr); 3547c6ae99SBarry Smith if (!aij) { 3647c6ae99SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)*J,"MatSeqAIJSetPreallocation_C",&aij);CHKERRQ(ierr); 3747c6ae99SBarry Smith } 3847c6ae99SBarry Smith if (aij) { 3947c6ae99SBarry Smith if (bs == 1) { 4047c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(*J,slice->d_nz,slice->d_nnz);CHKERRQ(ierr); 4147c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(*J,slice->d_nz,slice->d_nnz,slice->o_nz,slice->o_nnz);CHKERRQ(ierr); 4247c6ae99SBarry Smith } else if (!slice->d_nnz) { 43f77a5242SKarl Rupp ierr = MatSeqAIJSetPreallocation(*J,slice->d_nz*bs,NULL);CHKERRQ(ierr); 44f77a5242SKarl Rupp ierr = MatMPIAIJSetPreallocation(*J,slice->d_nz*bs,NULL,slice->o_nz*bs,NULL);CHKERRQ(ierr); 4547c6ae99SBarry Smith } else { 460c010503SBarry Smith /* The user has provided preallocation per block-row, convert it to per scalar-row respecting DMSlicedSetBlockFills() if applicable */ 47dcca6d9dSJed Brown ierr = PetscMalloc2(slice->n*bs,&sd_nnz,(!!slice->o_nnz)*slice->n*bs,&so_nnz);CHKERRQ(ierr); 4847c6ae99SBarry Smith for (i=0; i<slice->n*bs; i++) { 4947c6ae99SBarry Smith sd_nnz[i] = (slice->d_nnz[i/bs]-1) * (slice->ofill ? slice->ofill->i[i%bs+1]-slice->ofill->i[i%bs] : bs) 5047c6ae99SBarry Smith + (slice->dfill ? slice->dfill->i[i%bs+1]-slice->dfill->i[i%bs] : bs); 5147c6ae99SBarry Smith if (so_nnz) { 5247c6ae99SBarry Smith so_nnz[i] = slice->o_nnz[i/bs] * (slice->ofill ? slice->ofill->i[i%bs+1]-slice->ofill->i[i%bs] : bs); 5347c6ae99SBarry Smith } 5447c6ae99SBarry Smith } 5547c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(*J,slice->d_nz*bs,sd_nnz);CHKERRQ(ierr); 5647c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(*J,slice->d_nz*bs,sd_nnz,slice->o_nz*bs,so_nnz);CHKERRQ(ierr); 5747c6ae99SBarry Smith ierr = PetscFree2(sd_nnz,so_nnz);CHKERRQ(ierr); 5847c6ae99SBarry Smith } 5947c6ae99SBarry Smith } 6047c6ae99SBarry Smith 6147c6ae99SBarry Smith /* Set up the local to global map. For the scalar map, we have to translate to entry-wise indexing instead of block-wise. */ 6245b6f7e9SBarry Smith ierr = PetscMalloc1(slice->n+slice->Nghosts,&globals);CHKERRQ(ierr); 63f77a5242SKarl Rupp ierr = MatGetOwnershipRange(*J,&rstart,NULL);CHKERRQ(ierr); 6445b6f7e9SBarry Smith for (i=0; i<slice->n; i++) globals[i] = rstart/bs + i; 658865f1eaSKarl Rupp 6645b6f7e9SBarry Smith for (i=0; i<slice->Nghosts; i++) { 6745b6f7e9SBarry Smith globals[slice->n+i] = slice->ghosts[i]; 6847c6ae99SBarry Smith } 6945b6f7e9SBarry Smith ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_SELF,bs,slice->n+slice->Nghosts,globals,PETSC_OWN_POINTER,&lmap);CHKERRQ(ierr); 70784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(*J,lmap,lmap);CHKERRQ(ierr); 71fcfd50ebSBarry Smith ierr = ISLocalToGlobalMappingDestroy(&lmap);CHKERRQ(ierr); 72c6b011d8SStefano Zampini ierr = MatSetDM(*J,dm);CHKERRQ(ierr); 7347c6ae99SBarry Smith PetscFunctionReturn(0); 7447c6ae99SBarry Smith } 7547c6ae99SBarry Smith 7647c6ae99SBarry Smith /*@C 770c010503SBarry Smith DMSlicedSetGhosts - Sets the global indices of other processes elements that will 7847c6ae99SBarry Smith be ghosts on this process 7947c6ae99SBarry Smith 8047c6ae99SBarry Smith Not Collective 8147c6ae99SBarry Smith 8247c6ae99SBarry Smith Input Parameters: 830c010503SBarry Smith + slice - the DM object 8447c6ae99SBarry Smith . bs - block size 8547c6ae99SBarry Smith . nlocal - number of local (owned, non-ghost) blocks 8647c6ae99SBarry Smith . Nghosts - number of ghost blocks on this process 8747c6ae99SBarry Smith - ghosts - global indices of each ghost block 8847c6ae99SBarry Smith 8947c6ae99SBarry Smith Level: advanced 9047c6ae99SBarry Smith 91240edb01SJed Brown .seealso DMDestroy(), DMCreateGlobalVector() 9247c6ae99SBarry Smith 9347c6ae99SBarry Smith @*/ 947087cfbeSBarry Smith PetscErrorCode DMSlicedSetGhosts(DM dm,PetscInt bs,PetscInt nlocal,PetscInt Nghosts,const PetscInt ghosts[]) 9547c6ae99SBarry Smith { 9647c6ae99SBarry Smith PetscErrorCode ierr; 9747c6ae99SBarry Smith DM_Sliced *slice = (DM_Sliced*)dm->data; 9847c6ae99SBarry Smith 9947c6ae99SBarry Smith PetscFunctionBegin; 1000c010503SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 10147c6ae99SBarry Smith ierr = PetscFree(slice->ghosts);CHKERRQ(ierr); 102785e854fSJed Brown ierr = PetscMalloc1(Nghosts,&slice->ghosts);CHKERRQ(ierr); 103580bdb30SBarry Smith ierr = PetscArraycpy(slice->ghosts,ghosts,Nghosts);CHKERRQ(ierr); 10447c6ae99SBarry Smith slice->bs = bs; 10547c6ae99SBarry Smith slice->n = nlocal; 10647c6ae99SBarry Smith slice->Nghosts = Nghosts; 10747c6ae99SBarry Smith PetscFunctionReturn(0); 10847c6ae99SBarry Smith } 10947c6ae99SBarry Smith 11047c6ae99SBarry Smith /*@C 1110c010503SBarry Smith DMSlicedSetPreallocation - sets the matrix memory preallocation for matrices computed by DMSliced 11247c6ae99SBarry Smith 11347c6ae99SBarry Smith Not Collective 11447c6ae99SBarry Smith 11547c6ae99SBarry Smith Input Parameters: 1160c010503SBarry Smith + slice - the DM object 11747c6ae99SBarry Smith . d_nz - number of block nonzeros per block row in diagonal portion of local 11847c6ae99SBarry Smith submatrix (same for all local rows) 11947c6ae99SBarry Smith . d_nnz - array containing the number of block nonzeros in the various block rows 12047c6ae99SBarry Smith of the in diagonal portion of the local (possibly different for each block 121f77a5242SKarl Rupp row) or NULL. 12247c6ae99SBarry Smith . o_nz - number of block nonzeros per block row in the off-diagonal portion of local 12347c6ae99SBarry Smith submatrix (same for all local rows). 12447c6ae99SBarry Smith - o_nnz - array containing the number of nonzeros in the various block rows of the 12547c6ae99SBarry Smith off-diagonal portion of the local submatrix (possibly different for 126f77a5242SKarl Rupp each block row) or NULL. 12747c6ae99SBarry Smith 12847c6ae99SBarry Smith Notes: 12947c6ae99SBarry Smith See MatMPIBAIJSetPreallocation() for more details on preallocation. If a scalar matrix (AIJ) is 1300c010503SBarry Smith obtained with DMSlicedGetMatrix(), the correct preallocation will be set, respecting DMSlicedSetBlockFills(). 13147c6ae99SBarry Smith 13247c6ae99SBarry Smith Level: advanced 13347c6ae99SBarry Smith 134240edb01SJed Brown .seealso DMDestroy(), DMCreateGlobalVector(), MatMPIAIJSetPreallocation(), 1350c010503SBarry Smith MatMPIBAIJSetPreallocation(), DMSlicedGetMatrix(), DMSlicedSetBlockFills() 13647c6ae99SBarry Smith 13747c6ae99SBarry Smith @*/ 1387087cfbeSBarry Smith PetscErrorCode DMSlicedSetPreallocation(DM dm,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[]) 13947c6ae99SBarry Smith { 14047c6ae99SBarry Smith DM_Sliced *slice = (DM_Sliced*)dm->data; 14147c6ae99SBarry Smith 14247c6ae99SBarry Smith PetscFunctionBegin; 14347c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 14447c6ae99SBarry Smith slice->d_nz = d_nz; 14547c6ae99SBarry Smith slice->d_nnz = (PetscInt*)d_nnz; 14647c6ae99SBarry Smith slice->o_nz = o_nz; 14747c6ae99SBarry Smith slice->o_nnz = (PetscInt*)o_nnz; 14847c6ae99SBarry Smith PetscFunctionReturn(0); 14947c6ae99SBarry Smith } 15047c6ae99SBarry Smith 1510c010503SBarry Smith static PetscErrorCode DMSlicedSetBlockFills_Private(PetscInt bs,const PetscInt *fill,DMSlicedBlockFills **inf) 15247c6ae99SBarry Smith { 15347c6ae99SBarry Smith PetscErrorCode ierr; 15447c6ae99SBarry Smith PetscInt i,j,nz,*fi,*fj; 1550c010503SBarry Smith DMSlicedBlockFills *f; 15647c6ae99SBarry Smith 15747c6ae99SBarry Smith PetscFunctionBegin; 15847c6ae99SBarry Smith PetscValidPointer(inf,3); 1592a808120SBarry Smith if (*inf) {ierr = PetscFree3(*inf,(*inf)->i,(*inf)->j);CHKERRQ(ierr);} 16047c6ae99SBarry Smith if (!fill) PetscFunctionReturn(0); 16147c6ae99SBarry Smith for (i=0,nz=0; i<bs*bs; i++) if (fill[i]) nz++; 162dcca6d9dSJed Brown ierr = PetscMalloc3(1,&f,bs+1,&fi,nz,&fj);CHKERRQ(ierr); 16347c6ae99SBarry Smith f->bs = bs; 16447c6ae99SBarry Smith f->nz = nz; 16547c6ae99SBarry Smith f->i = fi; 16647c6ae99SBarry Smith f->j = fj; 16747c6ae99SBarry Smith for (i=0,nz=0; i<bs; i++) { 16847c6ae99SBarry Smith fi[i] = nz; 16947c6ae99SBarry Smith for (j=0; j<bs; j++) if (fill[i*bs+j]) fj[nz++] = j; 17047c6ae99SBarry Smith } 17147c6ae99SBarry Smith fi[i] = nz; 17247c6ae99SBarry Smith *inf = f; 17347c6ae99SBarry Smith PetscFunctionReturn(0); 17447c6ae99SBarry Smith } 17547c6ae99SBarry Smith 17647c6ae99SBarry Smith /*@C 1770c010503SBarry Smith DMSlicedSetBlockFills - Sets the fill pattern in each block for a multi-component problem 1780c010503SBarry Smith of the matrix returned by DMSlicedGetMatrix(). 17947c6ae99SBarry Smith 180d083f849SBarry Smith Logically Collective on dm 18147c6ae99SBarry Smith 18247c6ae99SBarry Smith Input Parameter: 1830c010503SBarry Smith + sliced - the DM object 184f77a5242SKarl Rupp . dfill - the fill pattern in the diagonal block (may be NULL, means use dense block) 18547c6ae99SBarry Smith - ofill - the fill pattern in the off-diagonal blocks 18647c6ae99SBarry Smith 18747c6ae99SBarry Smith Notes: 18847c6ae99SBarry Smith This only makes sense for multicomponent problems using scalar matrix formats (AIJ). 189aa219208SBarry Smith See DMDASetBlockFills() for example usage. 19047c6ae99SBarry Smith 19147c6ae99SBarry Smith Level: advanced 19247c6ae99SBarry Smith 193aa219208SBarry Smith .seealso DMSlicedGetMatrix(), DMDASetBlockFills() 19447c6ae99SBarry Smith @*/ 1957087cfbeSBarry Smith PetscErrorCode DMSlicedSetBlockFills(DM dm,const PetscInt *dfill,const PetscInt *ofill) 19647c6ae99SBarry Smith { 19747c6ae99SBarry Smith DM_Sliced *slice = (DM_Sliced*)dm->data; 19847c6ae99SBarry Smith PetscErrorCode ierr; 19947c6ae99SBarry Smith 20047c6ae99SBarry Smith PetscFunctionBegin; 2010c010503SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2020c010503SBarry Smith ierr = DMSlicedSetBlockFills_Private(slice->bs,dfill,&slice->dfill);CHKERRQ(ierr); 2030c010503SBarry Smith ierr = DMSlicedSetBlockFills_Private(slice->bs,ofill,&slice->ofill);CHKERRQ(ierr); 20447c6ae99SBarry Smith PetscFunctionReturn(0); 20547c6ae99SBarry Smith } 20647c6ae99SBarry Smith 2073efe6655SBarry Smith static PetscErrorCode DMDestroy_Sliced(DM dm) 20847c6ae99SBarry Smith { 20947c6ae99SBarry Smith PetscErrorCode ierr; 21047c6ae99SBarry Smith DM_Sliced *slice = (DM_Sliced*)dm->data; 21147c6ae99SBarry Smith 21247c6ae99SBarry Smith PetscFunctionBegin; 21347c6ae99SBarry Smith ierr = PetscFree(slice->ghosts);CHKERRQ(ierr); 2143c4aec1bSBarry Smith if (slice->dfill) {ierr = PetscFree3(slice->dfill,slice->dfill->i,slice->dfill->j);CHKERRQ(ierr);} 2153c4aec1bSBarry Smith if (slice->ofill) {ierr = PetscFree3(slice->ofill,slice->ofill->i,slice->ofill->j);CHKERRQ(ierr);} 216435a35e8SMatthew G Knepley /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */ 217435a35e8SMatthew G Knepley ierr = PetscFree(slice);CHKERRQ(ierr); 21847c6ae99SBarry Smith PetscFunctionReturn(0); 21947c6ae99SBarry Smith } 22047c6ae99SBarry Smith 2213efe6655SBarry Smith static PetscErrorCode DMCreateGlobalVector_Sliced(DM dm,Vec *gvec) 22247c6ae99SBarry Smith { 22347c6ae99SBarry Smith PetscErrorCode ierr; 22447c6ae99SBarry Smith DM_Sliced *slice = (DM_Sliced*)dm->data; 22547c6ae99SBarry Smith 22647c6ae99SBarry Smith PetscFunctionBegin; 22747c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 22847c6ae99SBarry Smith PetscValidPointer(gvec,2); 229*ea78f98cSLisandro Dalcin *gvec = NULL; 230ce94432eSBarry Smith ierr = VecCreateGhostBlock(PetscObjectComm((PetscObject)dm),slice->bs,slice->n*slice->bs,PETSC_DETERMINE,slice->Nghosts,slice->ghosts,gvec);CHKERRQ(ierr); 231c688c046SMatthew G Knepley ierr = VecSetDM(*gvec,dm);CHKERRQ(ierr); 23247c6ae99SBarry Smith PetscFunctionReturn(0); 23347c6ae99SBarry Smith } 23447c6ae99SBarry Smith 2353efe6655SBarry Smith static PetscErrorCode DMGlobalToLocalBegin_Sliced(DM da,Vec g,InsertMode mode,Vec l) 2363efe6655SBarry Smith { 2373efe6655SBarry Smith PetscErrorCode ierr; 2383efe6655SBarry Smith PetscBool flg; 2393efe6655SBarry Smith 2403efe6655SBarry Smith PetscFunctionBegin; 241a3314f2cSMatthew G Knepley ierr = VecGhostIsLocalForm(g,l,&flg);CHKERRQ(ierr); 242ce94432eSBarry Smith if (!flg) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONG,"Local vector is not local form of global vector"); 2433efe6655SBarry Smith ierr = VecGhostUpdateEnd(g,mode,SCATTER_FORWARD);CHKERRQ(ierr); 2443efe6655SBarry Smith ierr = VecGhostUpdateBegin(g,mode,SCATTER_FORWARD);CHKERRQ(ierr); 2453efe6655SBarry Smith PetscFunctionReturn(0); 2463efe6655SBarry Smith } 2473efe6655SBarry Smith 2483efe6655SBarry Smith static PetscErrorCode DMGlobalToLocalEnd_Sliced(DM da,Vec g,InsertMode mode,Vec l) 2493efe6655SBarry Smith { 2503efe6655SBarry Smith PetscErrorCode ierr; 2513efe6655SBarry Smith PetscBool flg; 2523efe6655SBarry Smith 2533efe6655SBarry Smith PetscFunctionBegin; 2543efe6655SBarry Smith ierr = VecGhostIsLocalForm(g,l,&flg);CHKERRQ(ierr); 255ce94432eSBarry Smith if (!flg) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONG,"Local vector is not local form of global vector"); 2563efe6655SBarry Smith ierr = VecGhostUpdateEnd(g,mode,SCATTER_FORWARD);CHKERRQ(ierr); 2573efe6655SBarry Smith PetscFunctionReturn(0); 2583efe6655SBarry Smith } 2593efe6655SBarry Smith 2603e0de288SBarry Smith /*MC 2613e0de288SBarry Smith DMSLICED = "sliced" - A DM object that is used to manage data for a general graph. Uses VecCreateGhost() ghosted vectors for storing the fields 2623e0de288SBarry Smith 2633e0de288SBarry Smith See DMCreateSliced() for details. 2643e0de288SBarry Smith 2653e0de288SBarry Smith Level: intermediate 2663e0de288SBarry Smith 2673e0de288SBarry Smith .seealso: DMType, DMCOMPOSITE, DMCreateSliced(), DMCreate() 2683e0de288SBarry Smith M*/ 2693e0de288SBarry Smith 2708cc058d9SJed Brown PETSC_EXTERN PetscErrorCode DMCreate_Sliced(DM p) 271a4121054SBarry Smith { 272a4121054SBarry Smith PetscErrorCode ierr; 273a4121054SBarry Smith DM_Sliced *slice; 274a4121054SBarry Smith 275a4121054SBarry Smith PetscFunctionBegin; 276b00a9115SJed Brown ierr = PetscNewLog(p,&slice);CHKERRQ(ierr); 277a4121054SBarry Smith p->data = slice; 278a4121054SBarry Smith 279a4121054SBarry Smith p->ops->createglobalvector = DMCreateGlobalVector_Sliced; 28025296bd5SBarry Smith p->ops->creatematrix = DMCreateMatrix_Sliced; 2813efe6655SBarry Smith p->ops->globaltolocalbegin = DMGlobalToLocalBegin_Sliced; 2823efe6655SBarry Smith p->ops->globaltolocalend = DMGlobalToLocalEnd_Sliced; 283a4121054SBarry Smith p->ops->destroy = DMDestroy_Sliced; 284a4121054SBarry Smith PetscFunctionReturn(0); 285a4121054SBarry Smith } 286a4121054SBarry Smith 28747c6ae99SBarry Smith /*@C 2880c010503SBarry Smith DMSlicedCreate - Creates a DM object, used to manage data for a unstructured problem 28947c6ae99SBarry Smith 290d083f849SBarry Smith Collective 2910c010503SBarry Smith 2920c010503SBarry Smith Input Parameter: 2933efe6655SBarry Smith + comm - the processors that will share the global vector 2943efe6655SBarry Smith . bs - the block size 2953efe6655SBarry Smith . nlocal - number of vector entries on this process 2963efe6655SBarry Smith . Nghosts - number of ghost points needed on this process 2973efe6655SBarry Smith . ghosts - global indices of all ghost points for this process 2983efe6655SBarry Smith . d_nnz - matrix preallocation information representing coupling within this process 2993efe6655SBarry Smith - o_nnz - matrix preallocation information representing coupling between this process and other processes 3000c010503SBarry Smith 3010c010503SBarry Smith Output Parameters: 3020c010503SBarry Smith . slice - the slice object 3030c010503SBarry Smith 3043efe6655SBarry Smith Notes: 3053efe6655SBarry Smith This DM does not support DMCreateLocalVector(), DMGlobalToLocalBegin(), and DMGlobalToLocalEnd() instead one directly uses 3063efe6655SBarry Smith VecGhostGetLocalForm() and VecGhostRestoreLocalForm() to access the local representation and VecGhostUpdateBegin() and VecGhostUpdateEnd() to update 3073efe6655SBarry Smith the ghost points. 3083efe6655SBarry Smith 3093efe6655SBarry Smith One can use DMGlobalToLocalBegin(), and DMGlobalToLocalEnd() instead of VecGhostUpdateBegin() and VecGhostUpdateEnd(). 3103efe6655SBarry Smith 3110c010503SBarry Smith Level: advanced 3120c010503SBarry Smith 3133e0de288SBarry Smith .seealso DMDestroy(), DMCreateGlobalVector(), DMSetType(), DMSLICED, DMSlicedSetGhosts(), DMSlicedSetPreallocation(), VecGhostUpdateBegin(), VecGhostUpdateEnd(), 3143e0de288SBarry Smith VecGhostGetLocalForm(), VecGhostRestoreLocalForm() 3150c010503SBarry Smith 3160c010503SBarry Smith @*/ 3173efe6655SBarry 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) 3180c010503SBarry Smith { 3190c010503SBarry Smith PetscErrorCode ierr; 3200c010503SBarry Smith 3210c010503SBarry Smith PetscFunctionBegin; 3220c010503SBarry Smith PetscValidPointer(dm,2); 323a4121054SBarry Smith ierr = DMCreate(comm,dm);CHKERRQ(ierr); 324a4121054SBarry Smith ierr = DMSetType(*dm,DMSLICED);CHKERRQ(ierr); 3253efe6655SBarry Smith ierr = DMSlicedSetGhosts(*dm,bs,nlocal,Nghosts,ghosts);CHKERRQ(ierr); 3263efe6655SBarry Smith if (d_nnz) { 3273efe6655SBarry Smith ierr = DMSlicedSetPreallocation(*dm,0, d_nnz,0,o_nnz);CHKERRQ(ierr); 3283efe6655SBarry Smith } 3290c010503SBarry Smith PetscFunctionReturn(0); 3300c010503SBarry Smith } 3310c010503SBarry Smith 332