1c6db04a5SJed Brown #include <petscdmsliced.h> /*I "petscdmsliced.h" I*/ 207475bc1SBarry Smith #include <petscmat.h> 307475bc1SBarry 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 1647c6ae99SBarry Smith #undef __FUNCT__ 17950540a4SJed Brown #define __FUNCT__ "DMCreateMatrix_Sliced" 18b412c318SBarry Smith PetscErrorCode DMCreateMatrix_Sliced(DM dm, Mat *J) 1947c6ae99SBarry Smith { 2047c6ae99SBarry Smith PetscErrorCode ierr; 2147c6ae99SBarry Smith PetscInt *globals,*sd_nnz,*so_nnz,rstart,bs,i; 2247c6ae99SBarry Smith ISLocalToGlobalMapping lmap,blmap; 23f77a5242SKarl Rupp void (*aij)(void) = NULL; 2447c6ae99SBarry Smith DM_Sliced *slice = (DM_Sliced*)dm->data; 2547c6ae99SBarry Smith 2647c6ae99SBarry Smith PetscFunctionBegin; 2747c6ae99SBarry Smith bs = slice->bs; 28ce94432eSBarry Smith ierr = MatCreate(PetscObjectComm((PetscObject)dm),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); 31b412c318SBarry Smith ierr = MatSetType(*J,dm->mattype);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) { 45f77a5242SKarl Rupp ierr = MatSeqAIJSetPreallocation(*J,slice->d_nz*bs,NULL);CHKERRQ(ierr); 46f77a5242SKarl Rupp ierr = MatMPIAIJSetPreallocation(*J,slice->d_nz*bs,NULL,slice->o_nz*bs,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 */ 49dcca6d9dSJed Brown ierr = PetscMalloc2(slice->n*bs,&sd_nnz,(!!slice->o_nnz)*slice->n*bs,&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. */ 64785e854fSJed Brown ierr = PetscMalloc1((slice->n+slice->Nghosts)*bs,&globals);CHKERRQ(ierr); 65f77a5242SKarl Rupp ierr = MatGetOwnershipRange(*J,&rstart,NULL);CHKERRQ(ierr); 668865f1eaSKarl Rupp for (i=0; i<slice->n*bs; i++) globals[i] = rstart + i; 678865f1eaSKarl Rupp 6847c6ae99SBarry Smith for (i=0; i<slice->Nghosts*bs; i++) { 6947c6ae99SBarry Smith globals[slice->n*bs+i] = slice->ghosts[i/bs]*bs + i%bs; 7047c6ae99SBarry Smith } 7147c6ae99SBarry Smith ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_SELF,(slice->n+slice->Nghosts)*bs,globals,PETSC_OWN_POINTER,&lmap);CHKERRQ(ierr); 7247c6ae99SBarry Smith ierr = ISLocalToGlobalMappingBlock(lmap,bs,&blmap);CHKERRQ(ierr); 73784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(*J,lmap,lmap);CHKERRQ(ierr); 74784ac674SJed Brown ierr = MatSetLocalToGlobalMappingBlock(*J,blmap,blmap);CHKERRQ(ierr); 75fcfd50ebSBarry Smith ierr = ISLocalToGlobalMappingDestroy(&lmap);CHKERRQ(ierr); 76fcfd50ebSBarry Smith ierr = ISLocalToGlobalMappingDestroy(&blmap);CHKERRQ(ierr); 7747c6ae99SBarry Smith PetscFunctionReturn(0); 7847c6ae99SBarry Smith } 7947c6ae99SBarry Smith 8047c6ae99SBarry Smith #undef __FUNCT__ 810c010503SBarry Smith #define __FUNCT__ "DMSlicedSetGhosts" 8247c6ae99SBarry Smith /*@C 830c010503SBarry Smith DMSlicedSetGhosts - Sets the global indices of other processes elements that will 8447c6ae99SBarry Smith be ghosts on this process 8547c6ae99SBarry Smith 8647c6ae99SBarry Smith Not Collective 8747c6ae99SBarry Smith 8847c6ae99SBarry Smith Input Parameters: 890c010503SBarry Smith + slice - the DM object 9047c6ae99SBarry Smith . bs - block size 9147c6ae99SBarry Smith . nlocal - number of local (owned, non-ghost) blocks 9247c6ae99SBarry Smith . Nghosts - number of ghost blocks on this process 9347c6ae99SBarry Smith - ghosts - global indices of each ghost block 9447c6ae99SBarry Smith 9547c6ae99SBarry Smith Level: advanced 9647c6ae99SBarry Smith 97240edb01SJed Brown .seealso DMDestroy(), DMCreateGlobalVector() 9847c6ae99SBarry Smith 9947c6ae99SBarry Smith @*/ 1007087cfbeSBarry Smith PetscErrorCode DMSlicedSetGhosts(DM dm,PetscInt bs,PetscInt nlocal,PetscInt Nghosts,const PetscInt ghosts[]) 10147c6ae99SBarry Smith { 10247c6ae99SBarry Smith PetscErrorCode ierr; 10347c6ae99SBarry Smith DM_Sliced *slice = (DM_Sliced*)dm->data; 10447c6ae99SBarry Smith 10547c6ae99SBarry Smith PetscFunctionBegin; 1060c010503SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 10747c6ae99SBarry Smith ierr = PetscFree(slice->ghosts);CHKERRQ(ierr); 108785e854fSJed Brown ierr = PetscMalloc1(Nghosts,&slice->ghosts);CHKERRQ(ierr); 10947c6ae99SBarry Smith ierr = PetscMemcpy(slice->ghosts,ghosts,Nghosts*sizeof(PetscInt));CHKERRQ(ierr); 11047c6ae99SBarry Smith slice->bs = bs; 11147c6ae99SBarry Smith slice->n = nlocal; 11247c6ae99SBarry Smith slice->Nghosts = Nghosts; 11347c6ae99SBarry Smith PetscFunctionReturn(0); 11447c6ae99SBarry Smith } 11547c6ae99SBarry Smith 11647c6ae99SBarry Smith #undef __FUNCT__ 1170c010503SBarry Smith #define __FUNCT__ "DMSlicedSetPreallocation" 11847c6ae99SBarry Smith /*@C 1190c010503SBarry Smith DMSlicedSetPreallocation - sets the matrix memory preallocation for matrices computed by DMSliced 12047c6ae99SBarry Smith 12147c6ae99SBarry Smith Not Collective 12247c6ae99SBarry Smith 12347c6ae99SBarry Smith Input Parameters: 1240c010503SBarry Smith + slice - the DM object 12547c6ae99SBarry Smith . d_nz - number of block nonzeros per block row in diagonal portion of local 12647c6ae99SBarry Smith submatrix (same for all local rows) 12747c6ae99SBarry Smith . d_nnz - array containing the number of block nonzeros in the various block rows 12847c6ae99SBarry Smith of the in diagonal portion of the local (possibly different for each block 129f77a5242SKarl Rupp row) or NULL. 13047c6ae99SBarry Smith . o_nz - number of block nonzeros per block row in the off-diagonal portion of local 13147c6ae99SBarry Smith submatrix (same for all local rows). 13247c6ae99SBarry Smith - o_nnz - array containing the number of nonzeros in the various block rows of the 13347c6ae99SBarry Smith off-diagonal portion of the local submatrix (possibly different for 134f77a5242SKarl Rupp each block row) or NULL. 13547c6ae99SBarry Smith 13647c6ae99SBarry Smith Notes: 13747c6ae99SBarry Smith See MatMPIBAIJSetPreallocation() for more details on preallocation. If a scalar matrix (AIJ) is 1380c010503SBarry Smith obtained with DMSlicedGetMatrix(), the correct preallocation will be set, respecting DMSlicedSetBlockFills(). 13947c6ae99SBarry Smith 14047c6ae99SBarry Smith Level: advanced 14147c6ae99SBarry Smith 142240edb01SJed Brown .seealso DMDestroy(), DMCreateGlobalVector(), MatMPIAIJSetPreallocation(), 1430c010503SBarry Smith MatMPIBAIJSetPreallocation(), DMSlicedGetMatrix(), DMSlicedSetBlockFills() 14447c6ae99SBarry Smith 14547c6ae99SBarry Smith @*/ 1467087cfbeSBarry Smith PetscErrorCode DMSlicedSetPreallocation(DM dm,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[]) 14747c6ae99SBarry Smith { 14847c6ae99SBarry Smith DM_Sliced *slice = (DM_Sliced*)dm->data; 14947c6ae99SBarry Smith 15047c6ae99SBarry Smith PetscFunctionBegin; 15147c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 15247c6ae99SBarry Smith slice->d_nz = d_nz; 15347c6ae99SBarry Smith slice->d_nnz = (PetscInt*)d_nnz; 15447c6ae99SBarry Smith slice->o_nz = o_nz; 15547c6ae99SBarry Smith slice->o_nnz = (PetscInt*)o_nnz; 15647c6ae99SBarry Smith PetscFunctionReturn(0); 15747c6ae99SBarry Smith } 15847c6ae99SBarry Smith 15947c6ae99SBarry Smith #undef __FUNCT__ 1600c010503SBarry Smith #define __FUNCT__ "DMSlicedSetBlockFills_Private" 1610c010503SBarry Smith static PetscErrorCode DMSlicedSetBlockFills_Private(PetscInt bs,const PetscInt *fill,DMSlicedBlockFills **inf) 16247c6ae99SBarry Smith { 16347c6ae99SBarry Smith PetscErrorCode ierr; 16447c6ae99SBarry Smith PetscInt i,j,nz,*fi,*fj; 1650c010503SBarry Smith DMSlicedBlockFills *f; 16647c6ae99SBarry Smith 16747c6ae99SBarry Smith PetscFunctionBegin; 16847c6ae99SBarry Smith PetscValidPointer(inf,3); 16947c6ae99SBarry Smith if (*inf) {ierr = PetscFree3((*inf)->i,(*inf)->j,*inf);CHKERRQ(ierr);} 17047c6ae99SBarry Smith if (!fill) PetscFunctionReturn(0); 17147c6ae99SBarry Smith for (i=0,nz=0; i<bs*bs; i++) if (fill[i]) nz++; 172dcca6d9dSJed Brown ierr = PetscMalloc3(1,&f,bs+1,&fi,nz,&fj);CHKERRQ(ierr); 17347c6ae99SBarry Smith f->bs = bs; 17447c6ae99SBarry Smith f->nz = nz; 17547c6ae99SBarry Smith f->i = fi; 17647c6ae99SBarry Smith f->j = fj; 17747c6ae99SBarry Smith for (i=0,nz=0; i<bs; i++) { 17847c6ae99SBarry Smith fi[i] = nz; 17947c6ae99SBarry Smith for (j=0; j<bs; j++) if (fill[i*bs+j]) fj[nz++] = j; 18047c6ae99SBarry Smith } 18147c6ae99SBarry Smith fi[i] = nz; 18247c6ae99SBarry Smith *inf = f; 18347c6ae99SBarry Smith PetscFunctionReturn(0); 18447c6ae99SBarry Smith } 18547c6ae99SBarry Smith 18647c6ae99SBarry Smith #undef __FUNCT__ 1870c010503SBarry Smith #define __FUNCT__ "DMSlicedSetBlockFills" 18847c6ae99SBarry Smith /*@C 1890c010503SBarry Smith DMSlicedSetBlockFills - Sets the fill pattern in each block for a multi-component problem 1900c010503SBarry Smith of the matrix returned by DMSlicedGetMatrix(). 19147c6ae99SBarry Smith 1920c010503SBarry Smith Logically Collective on DM 19347c6ae99SBarry Smith 19447c6ae99SBarry Smith Input Parameter: 1950c010503SBarry Smith + sliced - the DM object 196f77a5242SKarl Rupp . dfill - the fill pattern in the diagonal block (may be NULL, means use dense block) 19747c6ae99SBarry Smith - ofill - the fill pattern in the off-diagonal blocks 19847c6ae99SBarry Smith 19947c6ae99SBarry Smith Notes: 20047c6ae99SBarry Smith This only makes sense for multicomponent problems using scalar matrix formats (AIJ). 201aa219208SBarry Smith See DMDASetBlockFills() for example usage. 20247c6ae99SBarry Smith 20347c6ae99SBarry Smith Level: advanced 20447c6ae99SBarry Smith 205aa219208SBarry Smith .seealso DMSlicedGetMatrix(), DMDASetBlockFills() 20647c6ae99SBarry Smith @*/ 2077087cfbeSBarry Smith PetscErrorCode DMSlicedSetBlockFills(DM dm,const PetscInt *dfill,const PetscInt *ofill) 20847c6ae99SBarry Smith { 20947c6ae99SBarry Smith DM_Sliced *slice = (DM_Sliced*)dm->data; 21047c6ae99SBarry Smith PetscErrorCode ierr; 21147c6ae99SBarry Smith 21247c6ae99SBarry Smith PetscFunctionBegin; 2130c010503SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2140c010503SBarry Smith ierr = DMSlicedSetBlockFills_Private(slice->bs,dfill,&slice->dfill);CHKERRQ(ierr); 2150c010503SBarry Smith ierr = DMSlicedSetBlockFills_Private(slice->bs,ofill,&slice->ofill);CHKERRQ(ierr); 21647c6ae99SBarry Smith PetscFunctionReturn(0); 21747c6ae99SBarry Smith } 21847c6ae99SBarry Smith 21947c6ae99SBarry Smith #undef __FUNCT__ 2200c010503SBarry Smith #define __FUNCT__ "DMDestroy_Sliced" 2213efe6655SBarry Smith static PetscErrorCode DMDestroy_Sliced(DM dm) 22247c6ae99SBarry Smith { 22347c6ae99SBarry Smith PetscErrorCode ierr; 22447c6ae99SBarry Smith DM_Sliced *slice = (DM_Sliced*)dm->data; 22547c6ae99SBarry Smith 22647c6ae99SBarry Smith PetscFunctionBegin; 22747c6ae99SBarry Smith ierr = PetscFree(slice->ghosts);CHKERRQ(ierr); 2283c4aec1bSBarry Smith if (slice->dfill) {ierr = PetscFree3(slice->dfill,slice->dfill->i,slice->dfill->j);CHKERRQ(ierr);} 2293c4aec1bSBarry Smith if (slice->ofill) {ierr = PetscFree3(slice->ofill,slice->ofill->i,slice->ofill->j);CHKERRQ(ierr);} 230435a35e8SMatthew G Knepley /* This was originally freed in DMDestroy(), but that prevents reference counting of backend objects */ 231435a35e8SMatthew G Knepley ierr = PetscFree(slice);CHKERRQ(ierr); 23247c6ae99SBarry Smith PetscFunctionReturn(0); 23347c6ae99SBarry Smith } 23447c6ae99SBarry Smith 23547c6ae99SBarry Smith #undef __FUNCT__ 2360c010503SBarry Smith #define __FUNCT__ "DMCreateGlobalVector_Sliced" 2373efe6655SBarry Smith static PetscErrorCode DMCreateGlobalVector_Sliced(DM dm,Vec *gvec) 23847c6ae99SBarry Smith { 23947c6ae99SBarry Smith PetscErrorCode ierr; 24047c6ae99SBarry Smith DM_Sliced *slice = (DM_Sliced*)dm->data; 24147c6ae99SBarry Smith 24247c6ae99SBarry Smith PetscFunctionBegin; 24347c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 24447c6ae99SBarry Smith PetscValidPointer(gvec,2); 24547c6ae99SBarry Smith *gvec = 0; 246ce94432eSBarry Smith ierr = VecCreateGhostBlock(PetscObjectComm((PetscObject)dm),slice->bs,slice->n*slice->bs,PETSC_DETERMINE,slice->Nghosts,slice->ghosts,gvec);CHKERRQ(ierr); 247c688c046SMatthew G Knepley ierr = VecSetDM(*gvec,dm);CHKERRQ(ierr); 24847c6ae99SBarry Smith PetscFunctionReturn(0); 24947c6ae99SBarry Smith } 25047c6ae99SBarry Smith 2513efe6655SBarry Smith #undef __FUNCT__ 2523efe6655SBarry Smith #define __FUNCT__ "DMGlobalToLocalBegin_Sliced" 2533efe6655SBarry Smith static PetscErrorCode DMGlobalToLocalBegin_Sliced(DM da,Vec g,InsertMode mode,Vec l) 2543efe6655SBarry Smith { 2553efe6655SBarry Smith PetscErrorCode ierr; 2563efe6655SBarry Smith PetscBool flg; 2573efe6655SBarry Smith 2583efe6655SBarry Smith PetscFunctionBegin; 259a3314f2cSMatthew G Knepley ierr = VecGhostIsLocalForm(g,l,&flg);CHKERRQ(ierr); 260ce94432eSBarry Smith if (!flg) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONG,"Local vector is not local form of global vector"); 2613efe6655SBarry Smith ierr = VecGhostUpdateEnd(g,mode,SCATTER_FORWARD);CHKERRQ(ierr); 2623efe6655SBarry Smith ierr = VecGhostUpdateBegin(g,mode,SCATTER_FORWARD);CHKERRQ(ierr); 2633efe6655SBarry Smith PetscFunctionReturn(0); 2643efe6655SBarry Smith } 2653efe6655SBarry Smith 2663efe6655SBarry Smith #undef __FUNCT__ 2673efe6655SBarry Smith #define __FUNCT__ "DMGlobalToLocalEnd_Sliced" 2683efe6655SBarry Smith static PetscErrorCode DMGlobalToLocalEnd_Sliced(DM da,Vec g,InsertMode mode,Vec l) 2693efe6655SBarry Smith { 2703efe6655SBarry Smith PetscErrorCode ierr; 2713efe6655SBarry Smith PetscBool flg; 2723efe6655SBarry Smith 2733efe6655SBarry Smith PetscFunctionBegin; 2743efe6655SBarry Smith ierr = VecGhostIsLocalForm(g,l,&flg);CHKERRQ(ierr); 275ce94432eSBarry Smith if (!flg) SETERRQ(PetscObjectComm((PetscObject)da),PETSC_ERR_ARG_WRONG,"Local vector is not local form of global vector"); 2763efe6655SBarry Smith ierr = VecGhostUpdateEnd(g,mode,SCATTER_FORWARD);CHKERRQ(ierr); 2773efe6655SBarry Smith PetscFunctionReturn(0); 2783efe6655SBarry Smith } 2793efe6655SBarry Smith 2803e0de288SBarry Smith /*MC 2813e0de288SBarry Smith DMSLICED = "sliced" - A DM object that is used to manage data for a general graph. Uses VecCreateGhost() ghosted vectors for storing the fields 2823e0de288SBarry Smith 2833e0de288SBarry Smith See DMCreateSliced() for details. 2843e0de288SBarry Smith 2853e0de288SBarry Smith Level: intermediate 2863e0de288SBarry Smith 2873e0de288SBarry Smith .seealso: DMType, DMCOMPOSITE, DMCreateSliced(), DMCreate() 2883e0de288SBarry Smith M*/ 2893e0de288SBarry Smith 290a4121054SBarry Smith #undef __FUNCT__ 291a4121054SBarry Smith #define __FUNCT__ "DMCreate_Sliced" 2928cc058d9SJed Brown PETSC_EXTERN PetscErrorCode DMCreate_Sliced(DM p) 293a4121054SBarry Smith { 294a4121054SBarry Smith PetscErrorCode ierr; 295a4121054SBarry Smith DM_Sliced *slice; 296a4121054SBarry Smith 297a4121054SBarry Smith PetscFunctionBegin; 298*b00a9115SJed Brown ierr = PetscNewLog(p,&slice);CHKERRQ(ierr); 299a4121054SBarry Smith p->data = slice; 300a4121054SBarry Smith 301b549a33dSJed Brown ierr = PetscObjectChangeTypeName((PetscObject)p,DMSLICED);CHKERRQ(ierr); 3028865f1eaSKarl Rupp 303a4121054SBarry Smith p->ops->createglobalvector = DMCreateGlobalVector_Sliced; 30425296bd5SBarry Smith p->ops->creatematrix = DMCreateMatrix_Sliced; 3053efe6655SBarry Smith p->ops->globaltolocalbegin = DMGlobalToLocalBegin_Sliced; 3063efe6655SBarry Smith p->ops->globaltolocalend = DMGlobalToLocalEnd_Sliced; 307a4121054SBarry Smith p->ops->destroy = DMDestroy_Sliced; 308a4121054SBarry Smith PetscFunctionReturn(0); 309a4121054SBarry Smith } 310a4121054SBarry Smith 31147c6ae99SBarry Smith #undef __FUNCT__ 3120c010503SBarry Smith #define __FUNCT__ "DMSlicedCreate" 31347c6ae99SBarry Smith /*@C 3140c010503SBarry Smith DMSlicedCreate - Creates a DM object, used to manage data for a unstructured problem 31547c6ae99SBarry Smith 3160c010503SBarry Smith Collective on MPI_Comm 3170c010503SBarry Smith 3180c010503SBarry Smith Input Parameter: 3193efe6655SBarry Smith + comm - the processors that will share the global vector 3203efe6655SBarry Smith . bs - the block size 3213efe6655SBarry Smith . nlocal - number of vector entries on this process 3223efe6655SBarry Smith . Nghosts - number of ghost points needed on this process 3233efe6655SBarry Smith . ghosts - global indices of all ghost points for this process 3243efe6655SBarry Smith . d_nnz - matrix preallocation information representing coupling within this process 3253efe6655SBarry Smith - o_nnz - matrix preallocation information representing coupling between this process and other processes 3260c010503SBarry Smith 3270c010503SBarry Smith Output Parameters: 3280c010503SBarry Smith . slice - the slice object 3290c010503SBarry Smith 3303efe6655SBarry Smith Notes: 3313efe6655SBarry Smith This DM does not support DMCreateLocalVector(), DMGlobalToLocalBegin(), and DMGlobalToLocalEnd() instead one directly uses 3323efe6655SBarry Smith VecGhostGetLocalForm() and VecGhostRestoreLocalForm() to access the local representation and VecGhostUpdateBegin() and VecGhostUpdateEnd() to update 3333efe6655SBarry Smith the ghost points. 3343efe6655SBarry Smith 3353efe6655SBarry Smith One can use DMGlobalToLocalBegin(), and DMGlobalToLocalEnd() instead of VecGhostUpdateBegin() and VecGhostUpdateEnd(). 3363efe6655SBarry Smith 3370c010503SBarry Smith Level: advanced 3380c010503SBarry Smith 3393e0de288SBarry Smith .seealso DMDestroy(), DMCreateGlobalVector(), DMSetType(), DMSLICED, DMSlicedSetGhosts(), DMSlicedSetPreallocation(), VecGhostUpdateBegin(), VecGhostUpdateEnd(), 3403e0de288SBarry Smith VecGhostGetLocalForm(), VecGhostRestoreLocalForm() 3410c010503SBarry Smith 3420c010503SBarry Smith @*/ 3433efe6655SBarry 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) 3440c010503SBarry Smith { 3450c010503SBarry Smith PetscErrorCode ierr; 3460c010503SBarry Smith 3470c010503SBarry Smith PetscFunctionBegin; 3480c010503SBarry Smith PetscValidPointer(dm,2); 349a4121054SBarry Smith ierr = DMCreate(comm,dm);CHKERRQ(ierr); 350a4121054SBarry Smith ierr = DMSetType(*dm,DMSLICED);CHKERRQ(ierr); 3513efe6655SBarry Smith ierr = DMSlicedSetGhosts(*dm,bs,nlocal,Nghosts,ghosts);CHKERRQ(ierr); 3523efe6655SBarry Smith if (d_nnz) { 3533efe6655SBarry Smith ierr = DMSlicedSetPreallocation(*dm,0, d_nnz,0,o_nnz);CHKERRQ(ierr); 3543efe6655SBarry Smith } 3550c010503SBarry Smith PetscFunctionReturn(0); 3560c010503SBarry Smith } 3570c010503SBarry Smith 358