147c6ae99SBarry Smith #define PETSCDM_DLL 247c6ae99SBarry Smith 347c6ae99SBarry Smith #include "petscda.h" /*I "petscda.h" I*/ 447c6ae99SBarry Smith #include "petscmat.h" /*I "petscmat.h" I*/ 547c6ae99SBarry Smith #include "private/dmimpl.h" /*I "petscmat.h" I*/ 647c6ae99SBarry Smith 747c6ae99SBarry Smith /* CSR storage of the nonzero structure of a bs*bs matrix */ 847c6ae99SBarry Smith typedef struct { 947c6ae99SBarry Smith PetscInt bs,nz,*i,*j; 10*0c010503SBarry Smith } DMSlicedBlockFills; 1147c6ae99SBarry Smith 1247c6ae99SBarry Smith typedef struct { 1347c6ae99SBarry Smith Vec globalvector; 1447c6ae99SBarry Smith PetscInt bs,n,N,Nghosts,*ghosts; 1547c6ae99SBarry Smith PetscInt d_nz,o_nz,*d_nnz,*o_nnz; 16*0c010503SBarry Smith DMSlicedBlockFills *dfill,*ofill; 1747c6ae99SBarry Smith } DM_Sliced; 1847c6ae99SBarry Smith 1947c6ae99SBarry Smith #undef __FUNCT__ 20*0c010503SBarry Smith #define __FUNCT__ "DMGetMatrix_Sliced" 21*0c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMGetMatrix_Sliced(DM dm, const MatType mtype,Mat *J) 2247c6ae99SBarry Smith { 2347c6ae99SBarry Smith PetscErrorCode ierr; 2447c6ae99SBarry Smith PetscInt *globals,*sd_nnz,*so_nnz,rstart,bs,i; 2547c6ae99SBarry Smith ISLocalToGlobalMapping lmap,blmap; 2647c6ae99SBarry Smith void (*aij)(void) = PETSC_NULL; 2747c6ae99SBarry Smith DM_Sliced *slice = (DM_Sliced*)dm->data; 2847c6ae99SBarry Smith 2947c6ae99SBarry Smith PetscFunctionBegin; 3047c6ae99SBarry Smith bs = slice->bs; 3147c6ae99SBarry Smith ierr = MatCreate(((PetscObject)dm)->comm,J);CHKERRQ(ierr); 3247c6ae99SBarry Smith ierr = MatSetSizes(*J,slice->n*bs,slice->n*bs,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 3347c6ae99SBarry Smith ierr = MatSetType(*J,mtype);CHKERRQ(ierr); 3447c6ae99SBarry Smith ierr = MatSeqBAIJSetPreallocation(*J,bs,slice->d_nz,slice->d_nnz);CHKERRQ(ierr); 3547c6ae99SBarry Smith ierr = MatMPIBAIJSetPreallocation(*J,bs,slice->d_nz,slice->d_nnz,slice->o_nz,slice->o_nnz);CHKERRQ(ierr); 3647c6ae99SBarry Smith /* In general, we have to do extra work to preallocate for scalar (AIJ) matrices so we check whether it will do any 3747c6ae99SBarry Smith * good before going on with it. */ 3847c6ae99SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)*J,"MatMPIAIJSetPreallocation_C",&aij);CHKERRQ(ierr); 3947c6ae99SBarry Smith if (!aij) { 4047c6ae99SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)*J,"MatSeqAIJSetPreallocation_C",&aij);CHKERRQ(ierr); 4147c6ae99SBarry Smith } 4247c6ae99SBarry Smith if (aij) { 4347c6ae99SBarry Smith if (bs == 1) { 4447c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(*J,slice->d_nz,slice->d_nnz);CHKERRQ(ierr); 4547c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(*J,slice->d_nz,slice->d_nnz,slice->o_nz,slice->o_nnz);CHKERRQ(ierr); 4647c6ae99SBarry Smith } else if (!slice->d_nnz) { 4747c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(*J,slice->d_nz*bs,PETSC_NULL);CHKERRQ(ierr); 4847c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(*J,slice->d_nz*bs,PETSC_NULL,slice->o_nz*bs,PETSC_NULL);CHKERRQ(ierr); 4947c6ae99SBarry Smith } else { 50*0c010503SBarry Smith /* The user has provided preallocation per block-row, convert it to per scalar-row respecting DMSlicedSetBlockFills() if applicable */ 5147c6ae99SBarry Smith ierr = PetscMalloc2(slice->n*bs,PetscInt,&sd_nnz,(!!slice->o_nnz)*slice->n*bs,PetscInt,&so_nnz);CHKERRQ(ierr); 5247c6ae99SBarry Smith for (i=0; i<slice->n*bs; i++) { 5347c6ae99SBarry Smith sd_nnz[i] = (slice->d_nnz[i/bs]-1) * (slice->ofill ? slice->ofill->i[i%bs+1]-slice->ofill->i[i%bs] : bs) 5447c6ae99SBarry Smith + (slice->dfill ? slice->dfill->i[i%bs+1]-slice->dfill->i[i%bs] : bs); 5547c6ae99SBarry Smith if (so_nnz) { 5647c6ae99SBarry Smith so_nnz[i] = slice->o_nnz[i/bs] * (slice->ofill ? slice->ofill->i[i%bs+1]-slice->ofill->i[i%bs] : bs); 5747c6ae99SBarry Smith } 5847c6ae99SBarry Smith } 5947c6ae99SBarry Smith ierr = MatSeqAIJSetPreallocation(*J,slice->d_nz*bs,sd_nnz);CHKERRQ(ierr); 6047c6ae99SBarry Smith ierr = MatMPIAIJSetPreallocation(*J,slice->d_nz*bs,sd_nnz,slice->o_nz*bs,so_nnz);CHKERRQ(ierr); 6147c6ae99SBarry Smith ierr = PetscFree2(sd_nnz,so_nnz);CHKERRQ(ierr); 6247c6ae99SBarry Smith } 6347c6ae99SBarry Smith } 6447c6ae99SBarry Smith 6547c6ae99SBarry Smith ierr = MatSetBlockSize(*J,bs);CHKERRQ(ierr); 6647c6ae99SBarry Smith 6747c6ae99SBarry Smith /* Set up the local to global map. For the scalar map, we have to translate to entry-wise indexing instead of block-wise. */ 6847c6ae99SBarry Smith ierr = PetscMalloc((slice->n+slice->Nghosts)*bs*sizeof(PetscInt),&globals);CHKERRQ(ierr); 6947c6ae99SBarry Smith ierr = MatGetOwnershipRange(*J,&rstart,PETSC_NULL);CHKERRQ(ierr); 7047c6ae99SBarry Smith for (i=0; i<slice->n*bs; i++) { 7147c6ae99SBarry Smith globals[i] = rstart + i; 7247c6ae99SBarry Smith } 7347c6ae99SBarry Smith for (i=0; i<slice->Nghosts*bs; i++) { 7447c6ae99SBarry Smith globals[slice->n*bs+i] = slice->ghosts[i/bs]*bs + i%bs; 7547c6ae99SBarry Smith } 7647c6ae99SBarry Smith ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_SELF,(slice->n+slice->Nghosts)*bs,globals,PETSC_OWN_POINTER,&lmap);CHKERRQ(ierr); 7747c6ae99SBarry Smith ierr = ISLocalToGlobalMappingBlock(lmap,bs,&blmap);CHKERRQ(ierr); 7847c6ae99SBarry Smith ierr = MatSetLocalToGlobalMapping(*J,lmap);CHKERRQ(ierr); 7947c6ae99SBarry Smith ierr = MatSetLocalToGlobalMappingBlock(*J,blmap);CHKERRQ(ierr); 8047c6ae99SBarry Smith ierr = ISLocalToGlobalMappingDestroy(lmap);CHKERRQ(ierr); 8147c6ae99SBarry Smith ierr = ISLocalToGlobalMappingDestroy(blmap);CHKERRQ(ierr); 8247c6ae99SBarry Smith PetscFunctionReturn(0); 8347c6ae99SBarry Smith } 8447c6ae99SBarry Smith 8547c6ae99SBarry Smith #undef __FUNCT__ 86*0c010503SBarry Smith #define __FUNCT__ "DMSlicedSetGhosts" 8747c6ae99SBarry Smith /*@C 88*0c010503SBarry Smith DMSlicedSetGhosts - Sets the global indices of other processes elements that will 8947c6ae99SBarry Smith be ghosts on this process 9047c6ae99SBarry Smith 9147c6ae99SBarry Smith Not Collective 9247c6ae99SBarry Smith 9347c6ae99SBarry Smith Input Parameters: 94*0c010503SBarry Smith + slice - the DM object 9547c6ae99SBarry Smith . bs - block size 9647c6ae99SBarry Smith . nlocal - number of local (owned, non-ghost) blocks 9747c6ae99SBarry Smith . Nghosts - number of ghost blocks on this process 9847c6ae99SBarry Smith - ghosts - global indices of each ghost block 9947c6ae99SBarry Smith 10047c6ae99SBarry Smith Level: advanced 10147c6ae99SBarry Smith 102*0c010503SBarry Smith .seealso DMDestroy(), DMCreateGlobalVector(), DMSlicedGetGlobalIndices() 10347c6ae99SBarry Smith 10447c6ae99SBarry Smith @*/ 105*0c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMSlicedSetGhosts(DM dm,PetscInt bs,PetscInt nlocal,PetscInt Nghosts,const PetscInt ghosts[]) 10647c6ae99SBarry Smith { 10747c6ae99SBarry Smith PetscErrorCode ierr; 10847c6ae99SBarry Smith DM_Sliced *slice = (DM_Sliced*)dm->data; 10947c6ae99SBarry Smith 11047c6ae99SBarry Smith PetscFunctionBegin; 111*0c010503SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 11247c6ae99SBarry Smith ierr = PetscFree(slice->ghosts);CHKERRQ(ierr); 11347c6ae99SBarry Smith ierr = PetscMalloc(Nghosts*sizeof(PetscInt),&slice->ghosts);CHKERRQ(ierr); 11447c6ae99SBarry Smith ierr = PetscMemcpy(slice->ghosts,ghosts,Nghosts*sizeof(PetscInt));CHKERRQ(ierr); 11547c6ae99SBarry Smith slice->bs = bs; 11647c6ae99SBarry Smith slice->n = nlocal; 11747c6ae99SBarry Smith slice->Nghosts = Nghosts; 11847c6ae99SBarry Smith PetscFunctionReturn(0); 11947c6ae99SBarry Smith } 12047c6ae99SBarry Smith 12147c6ae99SBarry Smith #undef __FUNCT__ 122*0c010503SBarry Smith #define __FUNCT__ "DMSlicedSetPreallocation" 12347c6ae99SBarry Smith /*@C 124*0c010503SBarry Smith DMSlicedSetPreallocation - sets the matrix memory preallocation for matrices computed by DMSliced 12547c6ae99SBarry Smith 12647c6ae99SBarry Smith Not Collective 12747c6ae99SBarry Smith 12847c6ae99SBarry Smith Input Parameters: 129*0c010503SBarry Smith + slice - the DM object 13047c6ae99SBarry Smith . d_nz - number of block nonzeros per block row in diagonal portion of local 13147c6ae99SBarry Smith submatrix (same for all local rows) 13247c6ae99SBarry Smith . d_nnz - array containing the number of block nonzeros in the various block rows 13347c6ae99SBarry Smith of the in diagonal portion of the local (possibly different for each block 13447c6ae99SBarry Smith row) or PETSC_NULL. You must leave room for the diagonal entry even if it is zero. 13547c6ae99SBarry Smith . o_nz - number of block nonzeros per block row in the off-diagonal portion of local 13647c6ae99SBarry Smith submatrix (same for all local rows). 13747c6ae99SBarry Smith - o_nnz - array containing the number of nonzeros in the various block rows of the 13847c6ae99SBarry Smith off-diagonal portion of the local submatrix (possibly different for 13947c6ae99SBarry Smith each block row) or PETSC_NULL. 14047c6ae99SBarry Smith 14147c6ae99SBarry Smith Notes: 14247c6ae99SBarry Smith See MatMPIBAIJSetPreallocation() for more details on preallocation. If a scalar matrix (AIJ) is 143*0c010503SBarry Smith obtained with DMSlicedGetMatrix(), the correct preallocation will be set, respecting DMSlicedSetBlockFills(). 14447c6ae99SBarry Smith 14547c6ae99SBarry Smith Level: advanced 14647c6ae99SBarry Smith 147*0c010503SBarry Smith .seealso DMDestroy(), DMCreateGlobalVector(), DMSlicedGetGlobalIndices(), MatMPIAIJSetPreallocation(), 148*0c010503SBarry Smith MatMPIBAIJSetPreallocation(), DMSlicedGetMatrix(), DMSlicedSetBlockFills() 14947c6ae99SBarry Smith 15047c6ae99SBarry Smith @*/ 151*0c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMSlicedSetPreallocation(DM dm,PetscInt d_nz,const PetscInt d_nnz[],PetscInt o_nz,const PetscInt o_nnz[]) 15247c6ae99SBarry Smith { 15347c6ae99SBarry Smith DM_Sliced *slice = (DM_Sliced*)dm->data; 15447c6ae99SBarry Smith 15547c6ae99SBarry Smith PetscFunctionBegin; 15647c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 15747c6ae99SBarry Smith slice->d_nz = d_nz; 15847c6ae99SBarry Smith slice->d_nnz = (PetscInt*)d_nnz; 15947c6ae99SBarry Smith slice->o_nz = o_nz; 16047c6ae99SBarry Smith slice->o_nnz = (PetscInt*)o_nnz; 16147c6ae99SBarry Smith PetscFunctionReturn(0); 16247c6ae99SBarry Smith } 16347c6ae99SBarry Smith 16447c6ae99SBarry Smith #undef __FUNCT__ 165*0c010503SBarry Smith #define __FUNCT__ "DMSlicedSetBlockFills_Private" 166*0c010503SBarry Smith static PetscErrorCode DMSlicedSetBlockFills_Private(PetscInt bs,const PetscInt *fill,DMSlicedBlockFills **inf) 16747c6ae99SBarry Smith { 16847c6ae99SBarry Smith PetscErrorCode ierr; 16947c6ae99SBarry Smith PetscInt i,j,nz,*fi,*fj; 170*0c010503SBarry Smith DMSlicedBlockFills *f; 17147c6ae99SBarry Smith 17247c6ae99SBarry Smith PetscFunctionBegin; 17347c6ae99SBarry Smith PetscValidPointer(inf,3); 17447c6ae99SBarry Smith if (*inf) {ierr = PetscFree3((*inf)->i,(*inf)->j,*inf);CHKERRQ(ierr);} 17547c6ae99SBarry Smith if (!fill) PetscFunctionReturn(0); 17647c6ae99SBarry Smith for (i=0,nz=0; i<bs*bs; i++) if (fill[i]) nz++; 177*0c010503SBarry Smith ierr = PetscMalloc3(1,DMSlicedBlockFills,&f,bs+1,PetscInt,&fi,nz,PetscInt,&fj);CHKERRQ(ierr); 17847c6ae99SBarry Smith f->bs = bs; 17947c6ae99SBarry Smith f->nz = nz; 18047c6ae99SBarry Smith f->i = fi; 18147c6ae99SBarry Smith f->j = fj; 18247c6ae99SBarry Smith for (i=0,nz=0; i<bs; i++) { 18347c6ae99SBarry Smith fi[i] = nz; 18447c6ae99SBarry Smith for (j=0; j<bs; j++) if (fill[i*bs+j]) fj[nz++] = j; 18547c6ae99SBarry Smith } 18647c6ae99SBarry Smith fi[i] = nz; 18747c6ae99SBarry Smith *inf = f; 18847c6ae99SBarry Smith PetscFunctionReturn(0); 18947c6ae99SBarry Smith } 19047c6ae99SBarry Smith 19147c6ae99SBarry Smith #undef __FUNCT__ 192*0c010503SBarry Smith #define __FUNCT__ "DMSlicedSetBlockFills" 19347c6ae99SBarry Smith /*@C 194*0c010503SBarry Smith DMSlicedSetBlockFills - Sets the fill pattern in each block for a multi-component problem 195*0c010503SBarry Smith of the matrix returned by DMSlicedGetMatrix(). 19647c6ae99SBarry Smith 197*0c010503SBarry Smith Logically Collective on DM 19847c6ae99SBarry Smith 19947c6ae99SBarry Smith Input Parameter: 200*0c010503SBarry Smith + sliced - the DM object 20147c6ae99SBarry Smith . dfill - the fill pattern in the diagonal block (may be PETSC_NULL, means use dense block) 20247c6ae99SBarry Smith - ofill - the fill pattern in the off-diagonal blocks 20347c6ae99SBarry Smith 20447c6ae99SBarry Smith Notes: 20547c6ae99SBarry Smith This only makes sense for multicomponent problems using scalar matrix formats (AIJ). 20647c6ae99SBarry Smith See DASetBlockFills() for example usage. 20747c6ae99SBarry Smith 20847c6ae99SBarry Smith Level: advanced 20947c6ae99SBarry Smith 210*0c010503SBarry Smith .seealso DMSlicedGetMatrix(), DASetBlockFills() 21147c6ae99SBarry Smith @*/ 212*0c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMSlicedSetBlockFills(DM dm,const PetscInt *dfill,const PetscInt *ofill) 21347c6ae99SBarry Smith { 21447c6ae99SBarry Smith DM_Sliced *slice = (DM_Sliced*)dm->data; 21547c6ae99SBarry Smith PetscErrorCode ierr; 21647c6ae99SBarry Smith 21747c6ae99SBarry Smith PetscFunctionBegin; 218*0c010503SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 219*0c010503SBarry Smith ierr = DMSlicedSetBlockFills_Private(slice->bs,dfill,&slice->dfill);CHKERRQ(ierr); 220*0c010503SBarry Smith ierr = DMSlicedSetBlockFills_Private(slice->bs,ofill,&slice->ofill);CHKERRQ(ierr); 22147c6ae99SBarry Smith PetscFunctionReturn(0); 22247c6ae99SBarry Smith } 22347c6ae99SBarry Smith 22447c6ae99SBarry Smith extern PetscErrorCode DMDestroy_Private(DM,PetscBool *); 22547c6ae99SBarry Smith 22647c6ae99SBarry Smith #undef __FUNCT__ 227*0c010503SBarry Smith #define __FUNCT__ "DMDestroy_Sliced" 228*0c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMDestroy_Sliced(DM dm) 22947c6ae99SBarry Smith { 23047c6ae99SBarry Smith PetscErrorCode ierr; 23147c6ae99SBarry Smith PetscBool done; 23247c6ae99SBarry Smith DM_Sliced *slice = (DM_Sliced*)dm->data; 23347c6ae99SBarry Smith 23447c6ae99SBarry Smith PetscFunctionBegin; 23547c6ae99SBarry Smith ierr = DMDestroy_Private(dm,&done);CHKERRQ(ierr); 23647c6ae99SBarry Smith if (!done) PetscFunctionReturn(0); 23747c6ae99SBarry Smith 23847c6ae99SBarry Smith if (slice->globalvector) {ierr = VecDestroy(slice->globalvector);CHKERRQ(ierr);} 23947c6ae99SBarry Smith ierr = PetscFree(slice->ghosts);CHKERRQ(ierr); 24047c6ae99SBarry Smith if (slice->dfill) {ierr = PetscFree3(slice->dfill,slice->dfill->i,slice->dfill->j);CHKERRQ(ierr);} 24147c6ae99SBarry Smith if (slice->ofill) {ierr = PetscFree3(slice->ofill,slice->ofill->i,slice->ofill->j);CHKERRQ(ierr);} 24247c6ae99SBarry Smith ierr = PetscFree(dm->data);CHKERRQ(ierr); 24347c6ae99SBarry Smith ierr = PetscHeaderDestroy(dm);CHKERRQ(ierr); 24447c6ae99SBarry Smith PetscFunctionReturn(0); 24547c6ae99SBarry Smith } 24647c6ae99SBarry Smith 24747c6ae99SBarry Smith 248*0c010503SBarry Smith 24947c6ae99SBarry Smith #undef __FUNCT__ 250*0c010503SBarry Smith #define __FUNCT__ "DMCreateGlobalVector_Sliced" 251*0c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMCreateGlobalVector_Sliced(DM dm,Vec *gvec) 25247c6ae99SBarry Smith { 25347c6ae99SBarry Smith PetscErrorCode ierr; 25447c6ae99SBarry Smith PetscInt bs,cnt; 25547c6ae99SBarry Smith DM_Sliced *slice = (DM_Sliced*)dm->data; 25647c6ae99SBarry Smith 25747c6ae99SBarry Smith PetscFunctionBegin; 25847c6ae99SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 25947c6ae99SBarry Smith PetscValidPointer(gvec,2); 26047c6ae99SBarry Smith *gvec = 0; 26147c6ae99SBarry Smith if (slice->globalvector) { 26247c6ae99SBarry Smith ierr = PetscObjectGetReference((PetscObject)slice->globalvector,&cnt);CHKERRQ(ierr); 26347c6ae99SBarry Smith if (cnt == 1) { /* Nobody else has a reference so we can just reference it and give it away */ 26447c6ae99SBarry Smith *gvec = slice->globalvector; 26547c6ae99SBarry Smith ierr = PetscObjectReference((PetscObject)*gvec);CHKERRQ(ierr); 26647c6ae99SBarry Smith ierr = VecZeroEntries(*gvec);CHKERRQ(ierr); 26747c6ae99SBarry Smith } else { /* Someone else has a reference so we duplicate the global vector */ 26847c6ae99SBarry Smith ierr = VecDuplicate(slice->globalvector,gvec);CHKERRQ(ierr); 26947c6ae99SBarry Smith } 27047c6ae99SBarry Smith } else { 27147c6ae99SBarry Smith bs = slice->bs; 27247c6ae99SBarry Smith ierr = VecCreateGhostBlock(((PetscObject)dm)->comm,bs,slice->n*bs,PETSC_DETERMINE,slice->Nghosts,slice->ghosts,&slice->globalvector);CHKERRQ(ierr); 27347c6ae99SBarry Smith *gvec = slice->globalvector; 27447c6ae99SBarry Smith ierr = PetscObjectReference((PetscObject)*gvec);CHKERRQ(ierr); 27547c6ae99SBarry Smith } 27647c6ae99SBarry Smith PetscFunctionReturn(0); 27747c6ae99SBarry Smith } 27847c6ae99SBarry Smith 27947c6ae99SBarry Smith #undef __FUNCT__ 280*0c010503SBarry Smith #define __FUNCT__ "DMSlicedCreate" 28147c6ae99SBarry Smith /*@C 282*0c010503SBarry Smith DMSlicedCreate - Creates a DM object, used to manage data for a unstructured problem 28347c6ae99SBarry Smith 284*0c010503SBarry Smith Collective on MPI_Comm 285*0c010503SBarry Smith 286*0c010503SBarry Smith Input Parameter: 287*0c010503SBarry Smith . comm - the processors that will share the global vector 288*0c010503SBarry Smith 289*0c010503SBarry Smith Output Parameters: 290*0c010503SBarry Smith . slice - the slice object 291*0c010503SBarry Smith 292*0c010503SBarry Smith Level: advanced 293*0c010503SBarry Smith 294*0c010503SBarry Smith .seealso DMDestroy(), DMCreateGlobalVector(), DMSlicedGetGlobalIndices() 295*0c010503SBarry Smith 296*0c010503SBarry Smith @*/ 297*0c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMSlicedCreate(MPI_Comm comm,DM *dm) 298*0c010503SBarry Smith { 299*0c010503SBarry Smith PetscErrorCode ierr; 300*0c010503SBarry Smith DM p; 301*0c010503SBarry Smith DM_Sliced *slice; 302*0c010503SBarry Smith 303*0c010503SBarry Smith PetscFunctionBegin; 304*0c010503SBarry Smith PetscValidPointer(dm,2); 305*0c010503SBarry Smith #ifndef PETSC_USE_DYNAMIC_LIBRARIES 306*0c010503SBarry Smith ierr = DMInitializePackage(PETSC_NULL);CHKERRQ(ierr); 307*0c010503SBarry Smith #endif 308*0c010503SBarry Smith 309*0c010503SBarry Smith ierr = PetscHeaderCreate(p,_p_DM,struct _DMOps,DM_CLASSID,0,"DM",comm,DMDestroy,DMView);CHKERRQ(ierr); 310*0c010503SBarry Smith ierr = PetscNewLog(p,DM_Sliced,&slice);CHKERRQ(ierr); 311*0c010503SBarry Smith p->data = slice; 312*0c010503SBarry Smith 313*0c010503SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)p,"Sliced");CHKERRQ(ierr); 314*0c010503SBarry Smith p->ops->createglobalvector = DMCreateGlobalVector_Sliced; 315*0c010503SBarry Smith p->ops->getmatrix = DMGetMatrix_Sliced; 316*0c010503SBarry Smith p->ops->destroy = DMDestroy_Sliced; 317*0c010503SBarry Smith *dm = p; 318*0c010503SBarry Smith PetscFunctionReturn(0); 319*0c010503SBarry Smith } 320*0c010503SBarry Smith 321*0c010503SBarry Smith 322*0c010503SBarry Smith #undef __FUNCT__ 323*0c010503SBarry Smith #define __FUNCT__ "DMSlicedGetGlobalIndices" 324*0c010503SBarry Smith /*@C 325*0c010503SBarry Smith DMSlicedGetGlobalIndices - Gets the global indices for all the local entries 326*0c010503SBarry Smith 327*0c010503SBarry Smith Collective on DM 32847c6ae99SBarry Smith 32947c6ae99SBarry Smith Input Parameter: 33047c6ae99SBarry Smith . slice - the slice object 33147c6ae99SBarry Smith 33247c6ae99SBarry Smith Output Parameters: 33347c6ae99SBarry Smith . idx - the individual indices for each packed vector/array 33447c6ae99SBarry Smith 33547c6ae99SBarry Smith Level: advanced 33647c6ae99SBarry Smith 33747c6ae99SBarry Smith Notes: 33847c6ae99SBarry Smith The idx parameters should be freed by the calling routine with PetscFree() 33947c6ae99SBarry Smith 340*0c010503SBarry Smith .seealso DMSlicedDestroy(), DMSlicedCreateGlobalVector(), DMSlicedCreate() 34147c6ae99SBarry Smith 34247c6ae99SBarry Smith @*/ 343*0c010503SBarry Smith PetscErrorCode PETSCDM_DLLEXPORT DMSlicedGetGlobalIndices(DM dm,PetscInt *idx[]) 34447c6ae99SBarry Smith { 34547c6ae99SBarry Smith PetscFunctionReturn(0); 34647c6ae99SBarry Smith } 34747c6ae99SBarry Smith 34847c6ae99SBarry Smith 34947c6ae99SBarry Smith /* Explanation of the missing functions for DA-style handling of the local vector: 35047c6ae99SBarry Smith 351*0c010503SBarry Smith DMSlicedCreateLocalVector() 352*0c010503SBarry Smith DMSlicedGlobalToLocalBegin() 353*0c010503SBarry Smith DMSlicedGlobalToLocalEnd() 35447c6ae99SBarry Smith 355*0c010503SBarry Smith There is no way to get the global form from a local form, so DMSlicedCreateLocalVector() is a memory leak without 356*0c010503SBarry Smith external accounting for the global vector. Also, DMSliced intends the user to work with the VecGhost interface since the 35747c6ae99SBarry Smith ghosts are already ordered after the owned entries. Contrast this to a DA where the local vector has a special 35847c6ae99SBarry Smith ordering described by the structured grid, hence it cannot share memory with the global form. For this reason, users 359*0c010503SBarry Smith of DMSliced should work with the global vector and use 36047c6ae99SBarry Smith 36147c6ae99SBarry Smith VecGhostGetLocalForm(), VecGhostRestoreLocalForm() 36247c6ae99SBarry Smith VecGhostUpdateBegin(), VecGhostUpdateEnd() 36347c6ae99SBarry Smith 36447c6ae99SBarry Smith rather than the missing DA-style functions. This is conceptually simpler and offers better performance than is 36547c6ae99SBarry Smith possible with the DA-style interface. 36647c6ae99SBarry Smith */ 367