1032b8ab6SVijay Mahadevan #include <petsc-private/dmmbimpl.h> /*I "petscdm.h" I*/ 2032b8ab6SVijay Mahadevan #include <petsc-private/vecimpl.h> /*I "petscdm.h" I*/ 3032b8ab6SVijay Mahadevan 4032b8ab6SVijay Mahadevan #include <petscdmmoab.h> 5032b8ab6SVijay Mahadevan 6032b8ab6SVijay Mahadevan static PetscErrorCode DMMoab_Compute_NNZ_From_Connectivity(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscBool); 7032b8ab6SVijay Mahadevan 8032b8ab6SVijay Mahadevan #undef __FUNCT__ 9032b8ab6SVijay Mahadevan #define __FUNCT__ "DMCreateMatrix_Moab" 10032b8ab6SVijay Mahadevan PetscErrorCode DMCreateMatrix_Moab(DM dm, MatType mtype,Mat *J) 11032b8ab6SVijay Mahadevan { 12032b8ab6SVijay Mahadevan PetscErrorCode ierr; 13032b8ab6SVijay Mahadevan ISLocalToGlobalMapping ltog; 14032b8ab6SVijay Mahadevan MatType mltype; 15032b8ab6SVijay Mahadevan PetscInt i,nloc,count,dof,innz,ionz,nsize; 16032b8ab6SVijay Mahadevan DM_Moab *dmmoab = (DM_Moab*)dm->data; 17032b8ab6SVijay Mahadevan moab::Range *range=dmmoab->vowned; 18032b8ab6SVijay Mahadevan PetscInt *gindices,*nnz,*onz; 19032b8ab6SVijay Mahadevan char *tmp=0; 20032b8ab6SVijay Mahadevan moab::ErrorCode merr; 21032b8ab6SVijay Mahadevan 22032b8ab6SVijay Mahadevan PetscFunctionBegin; 23032b8ab6SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 24032b8ab6SVijay Mahadevan PetscValidPointer(J,3); 25032b8ab6SVijay Mahadevan nloc = dmmoab->vowned->size() * dmmoab->bs; 26032b8ab6SVijay Mahadevan 27*db66d124SVijay Mahadevan /* create the Matrix and set its type as specified by user */ 28032b8ab6SVijay Mahadevan ierr = MatCreate(dmmoab->pcomm->comm(), J);CHKERRQ(ierr); 29032b8ab6SVijay Mahadevan ierr = MatSetSizes(*J, nloc, nloc, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 30032b8ab6SVijay Mahadevan ierr = MatSetType(*J, mtype);CHKERRQ(ierr); 31032b8ab6SVijay Mahadevan ierr = MatSetFromOptions(*J);CHKERRQ(ierr); 32032b8ab6SVijay Mahadevan 33*db66d124SVijay Mahadevan /* next, need to allocate the non-zero arrays to enable pre-allocation */ 34032b8ab6SVijay Mahadevan ierr = MatGetType(*J, &mltype);CHKERRQ(ierr); /* in case user overrode the default type from command-line, re-check the type */ 35032b8ab6SVijay Mahadevan ierr = PetscStrstr(mltype, "baij", &tmp);CHKERRQ(ierr); 36032b8ab6SVijay Mahadevan nsize = (tmp ? dmmoab->nloc:dmmoab->nloc*dmmoab->bs); 37032b8ab6SVijay Mahadevan 38032b8ab6SVijay Mahadevan /* allocate the nnz, onz arrays based on block size and local nodes */ 39032b8ab6SVijay Mahadevan ierr = PetscMalloc(nsize*sizeof(PetscInt),&nnz);CHKERRQ(ierr); 40032b8ab6SVijay Mahadevan ierr = PetscMemzero(nnz,sizeof(PetscInt)*nsize);CHKERRQ(ierr); 41032b8ab6SVijay Mahadevan ierr = PetscMalloc(nsize*sizeof(PetscInt),&onz);CHKERRQ(ierr); 42032b8ab6SVijay Mahadevan ierr = PetscMemzero(onz,sizeof(PetscInt)*nsize);CHKERRQ(ierr); 43032b8ab6SVijay Mahadevan 44032b8ab6SVijay Mahadevan /* compute the nonzero pattern based on MOAB connectivity data for local elements */ 45032b8ab6SVijay Mahadevan ierr = DMMoab_Compute_NNZ_From_Connectivity(dm,&innz,nnz,&ionz,onz,(tmp?PETSC_TRUE:PETSC_FALSE));CHKERRQ(ierr); 46032b8ab6SVijay Mahadevan 47032b8ab6SVijay Mahadevan if (dmmoab->bs > 1 && tmp) { 48*db66d124SVijay Mahadevan /* Block matrix created, now set local to global mapping */ 49032b8ab6SVijay Mahadevan ierr = PetscMalloc(sizeof(PetscInt)*range->size()*dmmoab->bs, &gindices);CHKERRQ(ierr); 50032b8ab6SVijay Mahadevan moab::Range::iterator iter; 51032b8ab6SVijay Mahadevan for(iter = range->begin(),count=0; iter != range->end(); iter++,count+=dmmoab->bs) { 52032b8ab6SVijay Mahadevan merr = dmmoab->mbiface->tag_get_data(dmmoab->ltog_tag,&(*iter),1,&dof);MBERRNM(merr); 53032b8ab6SVijay Mahadevan for(i=0; i<dmmoab->bs; ++i) 54032b8ab6SVijay Mahadevan gindices[count+i] = (dof)*dmmoab->bs+i; 55032b8ab6SVijay Mahadevan } 56032b8ab6SVijay Mahadevan 57032b8ab6SVijay Mahadevan ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_SELF,range->size(),gindices,PETSC_COPY_VALUES,<og);CHKERRQ(ierr); 58032b8ab6SVijay Mahadevan ierr = MatSetLocalToGlobalMappingBlock(*J,ltog,ltog);CHKERRQ(ierr); 59032b8ab6SVijay Mahadevan 60*db66d124SVijay Mahadevan /* Clean up */ 61032b8ab6SVijay Mahadevan ierr = ISLocalToGlobalMappingDestroy(<og);CHKERRQ(ierr); 62032b8ab6SVijay Mahadevan ierr = PetscFree(gindices);CHKERRQ(ierr); 63032b8ab6SVijay Mahadevan 64032b8ab6SVijay Mahadevan } else { 65032b8ab6SVijay Mahadevan ierr = MatSetBlockSize(*J, dmmoab->bs);CHKERRQ(ierr); 66032b8ab6SVijay Mahadevan } 67032b8ab6SVijay Mahadevan 68032b8ab6SVijay Mahadevan /* set preallocation based on different supported Mat types */ 69032b8ab6SVijay Mahadevan ierr = MatSeqAIJSetPreallocation(*J, innz, nnz);CHKERRQ(ierr); 70032b8ab6SVijay Mahadevan ierr = MatMPIAIJSetPreallocation(*J, innz, nnz, ionz, onz);CHKERRQ(ierr); 71032b8ab6SVijay Mahadevan ierr = MatSeqBAIJSetPreallocation(*J, dmmoab->bs, innz, nnz);CHKERRQ(ierr); 72032b8ab6SVijay Mahadevan ierr = MatMPIBAIJSetPreallocation(*J, dmmoab->bs, innz, nnz, ionz, onz);CHKERRQ(ierr); 73032b8ab6SVijay Mahadevan 74032b8ab6SVijay Mahadevan ierr = MatSetDM(*J, dm);CHKERRQ(ierr); 75032b8ab6SVijay Mahadevan ierr = PetscFree(nnz);CHKERRQ(ierr); 76032b8ab6SVijay Mahadevan ierr = PetscFree(onz);CHKERRQ(ierr); 77032b8ab6SVijay Mahadevan PetscFunctionReturn(0); 78032b8ab6SVijay Mahadevan } 79032b8ab6SVijay Mahadevan 80032b8ab6SVijay Mahadevan 81032b8ab6SVijay Mahadevan #undef __FUNCT__ 82032b8ab6SVijay Mahadevan #define __FUNCT__ "DMMoab_Compute_NNZ_From_Connectivity" 83032b8ab6SVijay Mahadevan PetscErrorCode DMMoab_Compute_NNZ_From_Connectivity(DM dm,PetscInt* innz,PetscInt* nnz,PetscInt* ionz,PetscInt* onz,PetscBool isbaij) 84032b8ab6SVijay Mahadevan { 85032b8ab6SVijay Mahadevan PetscErrorCode ierr; 86032b8ab6SVijay Mahadevan PetscInt i,j,k,nloc,count,vpere,bs,nsize,nghost_found; 87032b8ab6SVijay Mahadevan DM_Moab *dmmoab = (DM_Moab*)dm->data; 88032b8ab6SVijay Mahadevan const moab::EntityHandle *connect; 89032b8ab6SVijay Mahadevan moab::Range *vowned=dmmoab->vowned,*elocal=dmmoab->elocal,*eghost=dmmoab->eghost; 90032b8ab6SVijay Mahadevan moab::Range::iterator iter; 91032b8ab6SVijay Mahadevan PetscInt *vertex_ids,firstvtx,dof,offset; 92032b8ab6SVijay Mahadevan moab::Tag id_tag; 93032b8ab6SVijay Mahadevan moab::ErrorCode merr; 94032b8ab6SVijay Mahadevan 95032b8ab6SVijay Mahadevan PetscFunctionBegin; 96032b8ab6SVijay Mahadevan bs = dmmoab->bs; 97032b8ab6SVijay Mahadevan nloc = dmmoab->nloc; 98032b8ab6SVijay Mahadevan nsize = (isbaij ? nloc:nloc*bs); 99032b8ab6SVijay Mahadevan 100032b8ab6SVijay Mahadevan ierr = DMMoabGetLocalToGlobalTag(dm,&id_tag);CHKERRQ(ierr); 101032b8ab6SVijay Mahadevan merr = dmmoab->mbiface->tag_iterate(id_tag,dmmoab->vlocal->begin(),dmmoab->vlocal->end(), 102032b8ab6SVijay Mahadevan count,reinterpret_cast<void*&>(vertex_ids));MBERRNM(merr); 103032b8ab6SVijay Mahadevan 104032b8ab6SVijay Mahadevan merr = dmmoab->mbiface->tag_get_data(dmmoab->ltog_tag,&(*vowned->begin()),1,&firstvtx);MBERRNM(merr); 105032b8ab6SVijay Mahadevan 106032b8ab6SVijay Mahadevan for(iter = elocal->begin(),i=0; iter != elocal->end(); iter++,i++) { 107032b8ab6SVijay Mahadevan /* Get connectivity information in canonical ordering for the local element */ 108032b8ab6SVijay Mahadevan merr = dmmoab->mbiface->get_connectivity((*iter),connect,vpere);MBERRNM(merr); 109032b8ab6SVijay Mahadevan 110032b8ab6SVijay Mahadevan nghost_found=0; 111032b8ab6SVijay Mahadevan /* loop over vertices and update the number of connectivity */ 112032b8ab6SVijay Mahadevan for (j=0;j<vpere;j++) { 113032b8ab6SVijay Mahadevan moab::Range::const_iterator giter = dmmoab->vghost->find(connect[j]); 114032b8ab6SVijay Mahadevan if (giter != dmmoab->vghost->end()) nghost_found++; 115032b8ab6SVijay Mahadevan } 116032b8ab6SVijay Mahadevan 117032b8ab6SVijay Mahadevan /* loop over vertices and update the number of connectivity */ 118032b8ab6SVijay Mahadevan for (j=0;j<vpere;j++) { 119032b8ab6SVijay Mahadevan merr = dmmoab->mbiface->tag_get_data(dmmoab->ltog_tag,&connect[j],1,&dof);MBERRNM(merr); 120032b8ab6SVijay Mahadevan 121032b8ab6SVijay Mahadevan moab::Range::const_iterator giter = dmmoab->vghost->find(connect[j]); 122032b8ab6SVijay Mahadevan if (giter != dmmoab->vghost->end()) continue; 123032b8ab6SVijay Mahadevan 124032b8ab6SVijay Mahadevan /* if block format, then all the block data are local only */ 125032b8ab6SVijay Mahadevan if (!isbaij) { 126032b8ab6SVijay Mahadevan offset=(dof-firstvtx)*bs; 127032b8ab6SVijay Mahadevan for (k=0;k<bs;k++) { 128032b8ab6SVijay Mahadevan nnz[offset+k]+=vpere-1; 129032b8ab6SVijay Mahadevan onz[offset+k]+=nghost_found; /* found a ghost non-owned node */ 130032b8ab6SVijay Mahadevan } 131032b8ab6SVijay Mahadevan } 132032b8ab6SVijay Mahadevan else { 133032b8ab6SVijay Mahadevan nnz[(dof-firstvtx)]+=vpere-1; 134032b8ab6SVijay Mahadevan onz[(dof-firstvtx)]+=nghost_found; /* found a ghost non-owned node */ 135032b8ab6SVijay Mahadevan } 136032b8ab6SVijay Mahadevan } 137032b8ab6SVijay Mahadevan } 138032b8ab6SVijay Mahadevan 139032b8ab6SVijay Mahadevan for(iter = eghost->begin(),i=0; iter != eghost->end(); iter++,i++) { 140032b8ab6SVijay Mahadevan /* Get connectivity information in canonical ordering for the local element */ 141032b8ab6SVijay Mahadevan merr = dmmoab->mbiface->get_connectivity((*iter),connect,vpere);MBERRNM(merr); 142032b8ab6SVijay Mahadevan 143032b8ab6SVijay Mahadevan nghost_found=0; 144032b8ab6SVijay Mahadevan /* loop over vertices and update the number of connectivity */ 145032b8ab6SVijay Mahadevan for (j=0;j<vpere;j++) { 146032b8ab6SVijay Mahadevan moab::Range::const_iterator giter = dmmoab->vghost->find(connect[j]); 147032b8ab6SVijay Mahadevan if (giter != dmmoab->vghost->end()) nghost_found++; 148032b8ab6SVijay Mahadevan } 149032b8ab6SVijay Mahadevan 150032b8ab6SVijay Mahadevan /* loop over vertices and update the number of connectivity */ 151032b8ab6SVijay Mahadevan for (j=0;j<vpere;j++) { 152032b8ab6SVijay Mahadevan moab::Range::const_iterator giter = dmmoab->vghost->find(connect[j]); 153032b8ab6SVijay Mahadevan if (giter != dmmoab->vghost->end()) continue; 154032b8ab6SVijay Mahadevan 155032b8ab6SVijay Mahadevan merr = dmmoab->mbiface->tag_get_data(dmmoab->ltog_tag,&connect[j],1,&dof);MBERRNM(merr); 156032b8ab6SVijay Mahadevan if (dof-firstvtx < 0) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_LIB,"Invalid local vertex ordering. Found ghost node with ID less than first vertex: [%i, %i]",dof,firstvtx); 157032b8ab6SVijay Mahadevan 158032b8ab6SVijay Mahadevan /* if block format, then all the block data are local only */ 159032b8ab6SVijay Mahadevan if (!isbaij) { 160032b8ab6SVijay Mahadevan offset=(dof-firstvtx)*bs; 161032b8ab6SVijay Mahadevan for (k=0;k<bs;k++) 162032b8ab6SVijay Mahadevan onz[offset+k]+=nghost_found; 163032b8ab6SVijay Mahadevan } 164032b8ab6SVijay Mahadevan else { 165032b8ab6SVijay Mahadevan onz[(dof-firstvtx)]+=nghost_found; 166032b8ab6SVijay Mahadevan } 167032b8ab6SVijay Mahadevan } 168032b8ab6SVijay Mahadevan } 169032b8ab6SVijay Mahadevan 170032b8ab6SVijay Mahadevan if (innz) *innz=0; 171032b8ab6SVijay Mahadevan if (ionz) *ionz=0; 172032b8ab6SVijay Mahadevan for (i=0;i<nsize;i++) { 173032b8ab6SVijay Mahadevan nnz[i]+=1; /* self count the node */ 174032b8ab6SVijay Mahadevan if (!isbaij) { 175032b8ab6SVijay Mahadevan nnz[i]*=bs; 176032b8ab6SVijay Mahadevan onz[i]*=bs; 177032b8ab6SVijay Mahadevan } 178032b8ab6SVijay Mahadevan if (innz && nnz[i]>*innz) *innz=nnz[i]; 179032b8ab6SVijay Mahadevan if (ionz && onz[i]>*ionz) *ionz=onz[i]; 180032b8ab6SVijay Mahadevan } 181032b8ab6SVijay Mahadevan PetscFunctionReturn(0); 182032b8ab6SVijay Mahadevan } 183032b8ab6SVijay Mahadevan 184032b8ab6SVijay Mahadevan 185