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 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__ 170c010503SBarry Smith #define __FUNCT__ "DMGetMatrix_Sliced" 187087cfbeSBarry Smith PetscErrorCode DMGetMatrix_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); 3047c6ae99SBarry Smith ierr = MatSetType(*J,mtype);CHKERRQ(ierr); 3147c6ae99SBarry Smith ierr = MatSeqBAIJSetPreallocation(*J,bs,slice->d_nz,slice->d_nnz);CHKERRQ(ierr); 3247c6ae99SBarry Smith ierr = MatMPIBAIJSetPreallocation(*J,bs,slice->d_nz,slice->d_nnz,slice->o_nz,slice->o_nnz);CHKERRQ(ierr); 3347c6ae99SBarry Smith /* In general, we have to do extra work to preallocate for scalar (AIJ) matrices so we check whether it will do any 3447c6ae99SBarry Smith * good before going on with it. */ 3547c6ae99SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)*J,"MatMPIAIJSetPreallocation_C",&aij);CHKERRQ(ierr); 3647c6ae99SBarry Smith if (!aij) { 3747c6ae99SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)*J,"MatSeqAIJSetPreallocation_C",&aij);CHKERRQ(ierr); 3847c6ae99SBarry Smith } 3947c6ae99SBarry Smith if (aij) { 4047c6ae99SBarry Smith if (bs == 1) { 4147c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(*J,slice->d_nz,slice->d_nnz);CHKERRQ(ierr); 4247c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(*J,slice->d_nz,slice->d_nnz,slice->o_nz,slice->o_nnz);CHKERRQ(ierr); 4347c6ae99SBarry Smith } else if (!slice->d_nnz) { 4447c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(*J,slice->d_nz*bs,PETSC_NULL);CHKERRQ(ierr); 4547c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(*J,slice->d_nz*bs,PETSC_NULL,slice->o_nz*bs,PETSC_NULL);CHKERRQ(ierr); 4647c6ae99SBarry Smith } else { 470c010503SBarry Smith /* The user has provided preallocation per block-row, convert it to per scalar-row respecting DMSlicedSetBlockFills() if applicable */ 4847c6ae99SBarry Smith ierr = PetscMalloc2(slice->n*bs,PetscInt,&sd_nnz,(!!slice->o_nnz)*slice->n*bs,PetscInt,&so_nnz);CHKERRQ(ierr); 4947c6ae99SBarry Smith for (i=0; i<slice->n*bs; i++) { 5047c6ae99SBarry Smith sd_nnz[i] = (slice->d_nnz[i/bs]-1) * (slice->ofill ? slice->ofill->i[i%bs+1]-slice->ofill->i[i%bs] : bs) 5147c6ae99SBarry Smith + (slice->dfill ? slice->dfill->i[i%bs+1]-slice->dfill->i[i%bs] : bs); 5247c6ae99SBarry Smith if (so_nnz) { 5347c6ae99SBarry Smith so_nnz[i] = slice->o_nnz[i/bs] * (slice->ofill ? slice->ofill->i[i%bs+1]-slice->ofill->i[i%bs] : bs); 5447c6ae99SBarry Smith } 5547c6ae99SBarry Smith } 5647c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(*J,slice->d_nz*bs,sd_nnz);CHKERRQ(ierr); 5747c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(*J,slice->d_nz*bs,sd_nnz,slice->o_nz*bs,so_nnz);CHKERRQ(ierr); 5847c6ae99SBarry Smith ierr = PetscFree2(sd_nnz,so_nnz);CHKERRQ(ierr); 5947c6ae99SBarry Smith } 6047c6ae99SBarry Smith } 6147c6ae99SBarry Smith 6247c6ae99SBarry Smith ierr = MatSetBlockSize(*J,bs);CHKERRQ(ierr); 6347c6ae99SBarry Smith 6447c6ae99SBarry Smith /* Set up the local to global map. For the scalar map, we have to translate to entry-wise indexing instead of block-wise. */ 6547c6ae99SBarry Smith ierr = PetscMalloc((slice->n+slice->Nghosts)*bs*sizeof(PetscInt),&globals);CHKERRQ(ierr); 6647c6ae99SBarry Smith ierr = MatGetOwnershipRange(*J,&rstart,PETSC_NULL);CHKERRQ(ierr); 6747c6ae99SBarry Smith for (i=0; i<slice->n*bs; i++) { 6847c6ae99SBarry Smith globals[i] = rstart + i; 6947c6ae99SBarry Smith } 7047c6ae99SBarry Smith for (i=0; i<slice->Nghosts*bs; i++) { 7147c6ae99SBarry Smith globals[slice->n*bs+i] = slice->ghosts[i/bs]*bs + i%bs; 7247c6ae99SBarry Smith } 7347c6ae99SBarry Smith ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_SELF,(slice->n+slice->Nghosts)*bs,globals,PETSC_OWN_POINTER,&lmap);CHKERRQ(ierr); 7447c6ae99SBarry Smith ierr = ISLocalToGlobalMappingBlock(lmap,bs,&blmap);CHKERRQ(ierr); 75784ac674SJed Brown ierr = MatSetLocalToGlobalMapping(*J,lmap,lmap);CHKERRQ(ierr); 76784ac674SJed Brown ierr = MatSetLocalToGlobalMappingBlock(*J,blmap,blmap);CHKERRQ(ierr); 77fcfd50ebSBarry Smith ierr = ISLocalToGlobalMappingDestroy(&lmap);CHKERRQ(ierr); 78fcfd50ebSBarry Smith ierr = ISLocalToGlobalMappingDestroy(&blmap);CHKERRQ(ierr); 7947c6ae99SBarry Smith PetscFunctionReturn(0); 8047c6ae99SBarry Smith } 8147c6ae99SBarry Smith 8247c6ae99SBarry Smith #undef __FUNCT__ 830c010503SBarry Smith #define __FUNCT__ "DMSlicedSetGhosts" 8447c6ae99SBarry Smith /*@C 850c010503SBarry Smith DMSlicedSetGhosts - Sets the global indices of other processes elements that will 8647c6ae99SBarry Smith be ghosts on this process 8747c6ae99SBarry Smith 8847c6ae99SBarry Smith Not Collective 8947c6ae99SBarry Smith 9047c6ae99SBarry Smith Input Parameters: 910c010503SBarry Smith + slice - the DM object 9247c6ae99SBarry Smith . bs - block size 9347c6ae99SBarry Smith . nlocal - number of local (owned, non-ghost) blocks 9447c6ae99SBarry Smith . Nghosts - number of ghost blocks on this process 9547c6ae99SBarry Smith - ghosts - global indices of each ghost block 9647c6ae99SBarry Smith 9747c6ae99SBarry Smith Level: advanced 9847c6ae99SBarry Smith 990c010503SBarry Smith .seealso DMDestroy(), DMCreateGlobalVector(), DMSlicedGetGlobalIndices() 10047c6ae99SBarry Smith 10147c6ae99SBarry Smith @*/ 1027087cfbeSBarry Smith PetscErrorCode DMSlicedSetGhosts(DM dm,PetscInt bs,PetscInt nlocal,PetscInt Nghosts,const PetscInt ghosts[]) 10347c6ae99SBarry Smith { 10447c6ae99SBarry Smith PetscErrorCode ierr; 10547c6ae99SBarry Smith DM_Sliced *slice = (DM_Sliced*)dm->data; 10647c6ae99SBarry Smith 10747c6ae99SBarry Smith PetscFunctionBegin; 1080c010503SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 10947c6ae99SBarry Smith ierr = PetscFree(slice->ghosts);CHKERRQ(ierr); 11047c6ae99SBarry Smith ierr = PetscMalloc(Nghosts*sizeof(PetscInt),&slice->ghosts);CHKERRQ(ierr); 11147c6ae99SBarry Smith ierr = PetscMemcpy(slice->ghosts,ghosts,Nghosts*sizeof(PetscInt));CHKERRQ(ierr); 11247c6ae99SBarry Smith slice->bs = bs; 11347c6ae99SBarry Smith slice->n = nlocal; 11447c6ae99SBarry Smith slice->Nghosts = Nghosts; 11547c6ae99SBarry Smith PetscFunctionReturn(0); 11647c6ae99SBarry Smith } 11747c6ae99SBarry Smith 11847c6ae99SBarry Smith #undef __FUNCT__ 1190c010503SBarry Smith #define __FUNCT__ "DMSlicedSetPreallocation" 12047c6ae99SBarry Smith /*@C 1210c010503SBarry Smith DMSlicedSetPreallocation - sets the matrix memory preallocation for matrices computed by DMSliced 12247c6ae99SBarry Smith 12347c6ae99SBarry Smith Not Collective 12447c6ae99SBarry Smith 12547c6ae99SBarry Smith Input Parameters: 1260c010503SBarry Smith + slice - the DM object 12747c6ae99SBarry Smith . d_nz - number of block nonzeros per block row in diagonal portion of local 12847c6ae99SBarry Smith submatrix (same for all local rows) 12947c6ae99SBarry Smith . d_nnz - array containing the number of block nonzeros in the various block rows 13047c6ae99SBarry Smith of the in diagonal portion of the local (possibly different for each block 13147c6ae99SBarry Smith row) or PETSC_NULL. You must leave room for the diagonal entry even if it is zero. 13247c6ae99SBarry Smith . o_nz - number of block nonzeros per block row in the off-diagonal portion of local 13347c6ae99SBarry Smith submatrix (same for all local rows). 13447c6ae99SBarry Smith - o_nnz - array containing the number of nonzeros in the various block rows of the 13547c6ae99SBarry Smith off-diagonal portion of the local submatrix (possibly different for 13647c6ae99SBarry Smith each block row) or PETSC_NULL. 13747c6ae99SBarry Smith 13847c6ae99SBarry Smith Notes: 13947c6ae99SBarry Smith See MatMPIBAIJSetPreallocation() for more details on preallocation. If a scalar matrix (AIJ) is 1400c010503SBarry Smith obtained with DMSlicedGetMatrix(), the correct preallocation will be set, respecting DMSlicedSetBlockFills(). 14147c6ae99SBarry Smith 14247c6ae99SBarry Smith Level: advanced 14347c6ae99SBarry Smith 1440c010503SBarry Smith .seealso DMDestroy(), DMCreateGlobalVector(), DMSlicedGetGlobalIndices(), MatMPIAIJSetPreallocation(), 1450c010503SBarry Smith MatMPIBAIJSetPreallocation(), DMSlicedGetMatrix(), DMSlicedSetBlockFills() 14647c6ae99SBarry Smith 14747c6ae99SBarry Smith @*/ 1487087cfbeSBarry Smith PetscErrorCode DMSlicedSetPreallocation(DM dm,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[]) 14947c6ae99SBarry Smith { 15047c6ae99SBarry Smith DM_Sliced *slice = (DM_Sliced*)dm->data; 15147c6ae99SBarry Smith 15247c6ae99SBarry Smith PetscFunctionBegin; 15347c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 15447c6ae99SBarry Smith slice->d_nz = d_nz; 15547c6ae99SBarry Smith slice->d_nnz = (PetscInt*)d_nnz; 15647c6ae99SBarry Smith slice->o_nz = o_nz; 15747c6ae99SBarry Smith slice->o_nnz = (PetscInt*)o_nnz; 15847c6ae99SBarry Smith PetscFunctionReturn(0); 15947c6ae99SBarry Smith } 16047c6ae99SBarry Smith 16147c6ae99SBarry Smith #undef __FUNCT__ 1620c010503SBarry Smith #define __FUNCT__ "DMSlicedSetBlockFills_Private" 1630c010503SBarry Smith static PetscErrorCode DMSlicedSetBlockFills_Private(PetscInt bs,const PetscInt *fill,DMSlicedBlockFills **inf) 16447c6ae99SBarry Smith { 16547c6ae99SBarry Smith PetscErrorCode ierr; 16647c6ae99SBarry Smith PetscInt i,j,nz,*fi,*fj; 1670c010503SBarry Smith DMSlicedBlockFills *f; 16847c6ae99SBarry Smith 16947c6ae99SBarry Smith PetscFunctionBegin; 17047c6ae99SBarry Smith PetscValidPointer(inf,3); 17147c6ae99SBarry Smith if (*inf) {ierr = PetscFree3((*inf)->i,(*inf)->j,*inf);CHKERRQ(ierr);} 17247c6ae99SBarry Smith if (!fill) PetscFunctionReturn(0); 17347c6ae99SBarry Smith for (i=0,nz=0; i<bs*bs; i++) if (fill[i]) nz++; 1740c010503SBarry Smith ierr = PetscMalloc3(1,DMSlicedBlockFills,&f,bs+1,PetscInt,&fi,nz,PetscInt,&fj);CHKERRQ(ierr); 17547c6ae99SBarry Smith f->bs = bs; 17647c6ae99SBarry Smith f->nz = nz; 17747c6ae99SBarry Smith f->i = fi; 17847c6ae99SBarry Smith f->j = fj; 17947c6ae99SBarry Smith for (i=0,nz=0; i<bs; i++) { 18047c6ae99SBarry Smith fi[i] = nz; 18147c6ae99SBarry Smith for (j=0; j<bs; j++) if (fill[i*bs+j]) fj[nz++] = j; 18247c6ae99SBarry Smith } 18347c6ae99SBarry Smith fi[i] = nz; 18447c6ae99SBarry Smith *inf = f; 18547c6ae99SBarry Smith PetscFunctionReturn(0); 18647c6ae99SBarry Smith } 18747c6ae99SBarry Smith 18847c6ae99SBarry Smith #undef __FUNCT__ 1890c010503SBarry Smith #define __FUNCT__ "DMSlicedSetBlockFills" 19047c6ae99SBarry Smith /*@C 1910c010503SBarry Smith DMSlicedSetBlockFills - Sets the fill pattern in each block for a multi-component problem 1920c010503SBarry Smith of the matrix returned by DMSlicedGetMatrix(). 19347c6ae99SBarry Smith 1940c010503SBarry Smith Logically Collective on DM 19547c6ae99SBarry Smith 19647c6ae99SBarry Smith Input Parameter: 1970c010503SBarry Smith + sliced - the DM object 19847c6ae99SBarry Smith . dfill - the fill pattern in the diagonal block (may be PETSC_NULL, means use dense block) 19947c6ae99SBarry Smith - ofill - the fill pattern in the off-diagonal blocks 20047c6ae99SBarry Smith 20147c6ae99SBarry Smith Notes: 20247c6ae99SBarry Smith This only makes sense for multicomponent problems using scalar matrix formats (AIJ). 203aa219208SBarry Smith See DMDASetBlockFills() for example usage. 20447c6ae99SBarry Smith 20547c6ae99SBarry Smith Level: advanced 20647c6ae99SBarry Smith 207aa219208SBarry Smith .seealso DMSlicedGetMatrix(), DMDASetBlockFills() 20847c6ae99SBarry Smith @*/ 2097087cfbeSBarry Smith PetscErrorCode DMSlicedSetBlockFills(DM dm,const PetscInt *dfill,const PetscInt *ofill) 21047c6ae99SBarry Smith { 21147c6ae99SBarry Smith DM_Sliced *slice = (DM_Sliced*)dm->data; 21247c6ae99SBarry Smith PetscErrorCode ierr; 21347c6ae99SBarry Smith 21447c6ae99SBarry Smith PetscFunctionBegin; 2150c010503SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2160c010503SBarry Smith ierr = DMSlicedSetBlockFills_Private(slice->bs,dfill,&slice->dfill);CHKERRQ(ierr); 2170c010503SBarry Smith ierr = DMSlicedSetBlockFills_Private(slice->bs,ofill,&slice->ofill);CHKERRQ(ierr); 21847c6ae99SBarry Smith PetscFunctionReturn(0); 21947c6ae99SBarry Smith } 22047c6ae99SBarry Smith 22147c6ae99SBarry Smith #undef __FUNCT__ 2220c010503SBarry Smith #define __FUNCT__ "DMDestroy_Sliced" 2237087cfbeSBarry Smith PetscErrorCode DMDestroy_Sliced(DM dm) 22447c6ae99SBarry Smith { 22547c6ae99SBarry Smith PetscErrorCode ierr; 22647c6ae99SBarry Smith DM_Sliced *slice = (DM_Sliced*)dm->data; 22747c6ae99SBarry Smith 22847c6ae99SBarry Smith PetscFunctionBegin; 22947c6ae99SBarry Smith ierr = PetscFree(slice->ghosts);CHKERRQ(ierr); 2303c4aec1bSBarry Smith if (slice->dfill) {ierr = PetscFree3(slice->dfill,slice->dfill->i,slice->dfill->j);CHKERRQ(ierr);} 2313c4aec1bSBarry Smith if (slice->ofill) {ierr = PetscFree3(slice->ofill,slice->ofill->i,slice->ofill->j);CHKERRQ(ierr);} 23247c6ae99SBarry Smith PetscFunctionReturn(0); 23347c6ae99SBarry Smith } 23447c6ae99SBarry Smith 23547c6ae99SBarry Smith #undef __FUNCT__ 2360c010503SBarry Smith #define __FUNCT__ "DMCreateGlobalVector_Sliced" 2377087cfbeSBarry Smith 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; 246*9ff76d2cSJed Brown ierr = VecCreateGhostBlock(((PetscObject)dm)->comm,slice->bs,slice->n*slice->bs,PETSC_DETERMINE,slice->Nghosts,slice->ghosts,gvec);CHKERRQ(ierr); 247b549a33dSJed Brown ierr = PetscObjectCompose((PetscObject)*gvec,"DM",(PetscObject)dm);CHKERRQ(ierr); 24847c6ae99SBarry Smith PetscFunctionReturn(0); 24947c6ae99SBarry Smith } 25047c6ae99SBarry Smith 251a4121054SBarry Smith EXTERN_C_BEGIN 252a4121054SBarry Smith #undef __FUNCT__ 253a4121054SBarry Smith #define __FUNCT__ "DMCreate_Sliced" 2547087cfbeSBarry Smith PetscErrorCode DMCreate_Sliced(DM p) 255a4121054SBarry Smith { 256a4121054SBarry Smith PetscErrorCode ierr; 257a4121054SBarry Smith DM_Sliced *slice; 258a4121054SBarry Smith 259a4121054SBarry Smith PetscFunctionBegin; 260a4121054SBarry Smith ierr = PetscNewLog(p,DM_Sliced,&slice);CHKERRQ(ierr); 261a4121054SBarry Smith p->data = slice; 262a4121054SBarry Smith 263b549a33dSJed Brown ierr = PetscObjectChangeTypeName((PetscObject)p,DMSLICED);CHKERRQ(ierr); 264a4121054SBarry Smith p->ops->createglobalvector = DMCreateGlobalVector_Sliced; 265a4121054SBarry Smith p->ops->getmatrix = DMGetMatrix_Sliced; 266a4121054SBarry Smith p->ops->destroy = DMDestroy_Sliced; 267a4121054SBarry Smith PetscFunctionReturn(0); 268a4121054SBarry Smith } 269a4121054SBarry Smith EXTERN_C_END 270a4121054SBarry Smith 27147c6ae99SBarry Smith #undef __FUNCT__ 2720c010503SBarry Smith #define __FUNCT__ "DMSlicedCreate" 27347c6ae99SBarry Smith /*@C 2740c010503SBarry Smith DMSlicedCreate - Creates a DM object, used to manage data for a unstructured problem 27547c6ae99SBarry Smith 2760c010503SBarry Smith Collective on MPI_Comm 2770c010503SBarry Smith 2780c010503SBarry Smith Input Parameter: 2790c010503SBarry Smith . comm - the processors that will share the global vector 2800c010503SBarry Smith 2810c010503SBarry Smith Output Parameters: 2820c010503SBarry Smith . slice - the slice object 2830c010503SBarry Smith 2840c010503SBarry Smith Level: advanced 2850c010503SBarry Smith 2860c010503SBarry Smith .seealso DMDestroy(), DMCreateGlobalVector(), DMSlicedGetGlobalIndices() 2870c010503SBarry Smith 2880c010503SBarry Smith @*/ 2897087cfbeSBarry Smith PetscErrorCode DMSlicedCreate(MPI_Comm comm,DM *dm) 2900c010503SBarry Smith { 2910c010503SBarry Smith PetscErrorCode ierr; 2920c010503SBarry Smith 2930c010503SBarry Smith PetscFunctionBegin; 2940c010503SBarry Smith PetscValidPointer(dm,2); 295a4121054SBarry Smith ierr = DMCreate(comm,dm);CHKERRQ(ierr); 296a4121054SBarry Smith ierr = DMSetType(*dm,DMSLICED);CHKERRQ(ierr); 2970c010503SBarry Smith PetscFunctionReturn(0); 2980c010503SBarry Smith } 2990c010503SBarry Smith 3000c010503SBarry Smith 3010c010503SBarry Smith #undef __FUNCT__ 3020c010503SBarry Smith #define __FUNCT__ "DMSlicedGetGlobalIndices" 3030c010503SBarry Smith /*@C 3040c010503SBarry Smith DMSlicedGetGlobalIndices - Gets the global indices for all the local entries 3050c010503SBarry Smith 3060c010503SBarry Smith Collective on DM 30747c6ae99SBarry Smith 30847c6ae99SBarry Smith Input Parameter: 30947c6ae99SBarry Smith . slice - the slice object 31047c6ae99SBarry Smith 31147c6ae99SBarry Smith Output Parameters: 31247c6ae99SBarry Smith . idx - the individual indices for each packed vector/array 31347c6ae99SBarry Smith 31447c6ae99SBarry Smith Level: advanced 31547c6ae99SBarry Smith 31647c6ae99SBarry Smith Notes: 31747c6ae99SBarry Smith The idx parameters should be freed by the calling routine with PetscFree() 31847c6ae99SBarry Smith 319b549a33dSJed Brown .seealso DMSlicedDestroy(), DMCreateGlobalVector(), DMSlicedCreate() 32047c6ae99SBarry Smith 32147c6ae99SBarry Smith @*/ 3227087cfbeSBarry Smith PetscErrorCode DMSlicedGetGlobalIndices(DM dm,PetscInt *idx[]) 32347c6ae99SBarry Smith { 32447c6ae99SBarry Smith PetscFunctionReturn(0); 32547c6ae99SBarry Smith } 32647c6ae99SBarry Smith 32747c6ae99SBarry Smith 328aa219208SBarry Smith /* Explanation of the missing functions for DMDA-style handling of the local vector: 32947c6ae99SBarry Smith 3300c010503SBarry Smith DMSlicedCreateLocalVector() 3310c010503SBarry Smith DMSlicedGlobalToLocalBegin() 3320c010503SBarry Smith DMSlicedGlobalToLocalEnd() 33347c6ae99SBarry Smith 3340c010503SBarry Smith There is no way to get the global form from a local form, so DMSlicedCreateLocalVector() is a memory leak without 3350c010503SBarry Smith external accounting for the global vector. Also, DMSliced intends the user to work with the VecGhost interface since the 336aa219208SBarry Smith ghosts are already ordered after the owned entries. Contrast this to a DMDA where the local vector has a special 33747c6ae99SBarry Smith ordering described by the structured grid, hence it cannot share memory with the global form. For this reason, users 3380c010503SBarry Smith of DMSliced should work with the global vector and use 33947c6ae99SBarry Smith 34047c6ae99SBarry Smith VecGhostGetLocalForm(), VecGhostRestoreLocalForm() 34147c6ae99SBarry Smith VecGhostUpdateBegin(), VecGhostUpdateEnd() 34247c6ae99SBarry Smith 343aa219208SBarry Smith rather than the missing DMDA-style functions. This is conceptually simpler and offers better performance than is 344aa219208SBarry Smith possible with the DMDA-style interface. 34547c6ae99SBarry Smith */ 346