1*032b8ab6SVijay Mahadevan #include <petsc-private/dmmbimpl.h> /*I "petscdm.h" I*/ 2*032b8ab6SVijay Mahadevan #include <petsc-private/vecimpl.h> /*I "petscdm.h" I*/ 3*032b8ab6SVijay Mahadevan 4*032b8ab6SVijay Mahadevan #include <petscdmmoab.h> 5*032b8ab6SVijay Mahadevan #include <MBTagConventions.hpp> 6*032b8ab6SVijay Mahadevan #include <sstream> 7*032b8ab6SVijay Mahadevan 8*032b8ab6SVijay Mahadevan static PetscErrorCode DMMoab_Compute_NNZ_From_Connectivity(DM,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscBool); 9*032b8ab6SVijay Mahadevan 10*032b8ab6SVijay Mahadevan #undef __FUNCT__ 11*032b8ab6SVijay Mahadevan #define __FUNCT__ "DMCreateMatrix_Moab" 12*032b8ab6SVijay Mahadevan PetscErrorCode DMCreateMatrix_Moab(DM dm, MatType mtype,Mat *J) 13*032b8ab6SVijay Mahadevan { 14*032b8ab6SVijay Mahadevan PetscErrorCode ierr; 15*032b8ab6SVijay Mahadevan ISLocalToGlobalMapping ltog; 16*032b8ab6SVijay Mahadevan MatType mltype; 17*032b8ab6SVijay Mahadevan PetscInt i,nloc,count,dof,innz,ionz,nsize; 18*032b8ab6SVijay Mahadevan DM_Moab *dmmoab = (DM_Moab*)dm->data; 19*032b8ab6SVijay Mahadevan moab::Range *range=dmmoab->vowned; 20*032b8ab6SVijay Mahadevan PetscInt *gindices,*nnz,*onz; 21*032b8ab6SVijay Mahadevan char *tmp=0; 22*032b8ab6SVijay Mahadevan moab::ErrorCode merr; 23*032b8ab6SVijay Mahadevan 24*032b8ab6SVijay Mahadevan PetscFunctionBegin; 25*032b8ab6SVijay Mahadevan PetscValidHeaderSpecific(dm,DM_CLASSID,1); 26*032b8ab6SVijay Mahadevan PetscValidPointer(J,3); 27*032b8ab6SVijay Mahadevan nloc = dmmoab->vowned->size() * dmmoab->bs; 28*032b8ab6SVijay Mahadevan // ierr = DMMoabCreateMatrix(dm,PETSC_NULL,dmmoab->bs,dmmoab->vowned,PETSC_TRUE,PETSC_TRUE,J);CHKERRQ(ierr); 29*032b8ab6SVijay Mahadevan 30*032b8ab6SVijay Mahadevan ierr = MatCreate(dmmoab->pcomm->comm(), J);CHKERRQ(ierr); 31*032b8ab6SVijay Mahadevan ierr = MatSetSizes(*J, nloc, nloc, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 32*032b8ab6SVijay Mahadevan ierr = MatSetType(*J, mtype);CHKERRQ(ierr); 33*032b8ab6SVijay Mahadevan ierr = MatSetFromOptions(*J);CHKERRQ(ierr); 34*032b8ab6SVijay Mahadevan 35*032b8ab6SVijay Mahadevan /* allocate the non-zeros properly */ 36*032b8ab6SVijay Mahadevan ierr = MatGetType(*J, &mltype);CHKERRQ(ierr); /* in case user overrode the default type from command-line, re-check the type */ 37*032b8ab6SVijay Mahadevan ierr = PetscStrstr(mltype, "baij", &tmp);CHKERRQ(ierr); 38*032b8ab6SVijay Mahadevan nsize = (tmp ? dmmoab->nloc:dmmoab->nloc*dmmoab->bs); 39*032b8ab6SVijay Mahadevan 40*032b8ab6SVijay Mahadevan /* allocate the nnz, onz arrays based on block size and local nodes */ 41*032b8ab6SVijay Mahadevan ierr = PetscMalloc(nsize*sizeof(PetscInt),&nnz);CHKERRQ(ierr); 42*032b8ab6SVijay Mahadevan ierr = PetscMemzero(nnz,sizeof(PetscInt)*nsize);CHKERRQ(ierr); 43*032b8ab6SVijay Mahadevan ierr = PetscMalloc(nsize*sizeof(PetscInt),&onz);CHKERRQ(ierr); 44*032b8ab6SVijay Mahadevan ierr = PetscMemzero(onz,sizeof(PetscInt)*nsize);CHKERRQ(ierr); 45*032b8ab6SVijay Mahadevan 46*032b8ab6SVijay Mahadevan /* compute the nonzero pattern based on MOAB connectivity data for local elements */ 47*032b8ab6SVijay Mahadevan ierr = DMMoab_Compute_NNZ_From_Connectivity(dm,&innz,nnz,&ionz,onz,(tmp?PETSC_TRUE:PETSC_FALSE));CHKERRQ(ierr); 48*032b8ab6SVijay Mahadevan 49*032b8ab6SVijay Mahadevan if (dmmoab->bs > 1 && tmp) { 50*032b8ab6SVijay Mahadevan // Block matrix created, now set local to global mapping: 51*032b8ab6SVijay Mahadevan ierr = PetscMalloc(sizeof(PetscInt)*range->size()*dmmoab->bs, &gindices);CHKERRQ(ierr); 52*032b8ab6SVijay Mahadevan moab::Range::iterator iter; 53*032b8ab6SVijay Mahadevan for(iter = range->begin(),count=0; iter != range->end(); iter++,count+=dmmoab->bs) { 54*032b8ab6SVijay Mahadevan merr = dmmoab->mbiface->tag_get_data(dmmoab->ltog_tag,&(*iter),1,&dof);MBERRNM(merr); 55*032b8ab6SVijay Mahadevan for(i=0; i<dmmoab->bs; ++i) 56*032b8ab6SVijay Mahadevan gindices[count+i] = (dof)*dmmoab->bs+i; 57*032b8ab6SVijay Mahadevan } 58*032b8ab6SVijay Mahadevan 59*032b8ab6SVijay Mahadevan ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_SELF,range->size(),gindices,PETSC_COPY_VALUES,<og);CHKERRQ(ierr); 60*032b8ab6SVijay Mahadevan ierr = MatSetLocalToGlobalMappingBlock(*J,ltog,ltog);CHKERRQ(ierr); 61*032b8ab6SVijay Mahadevan 62*032b8ab6SVijay Mahadevan // Clean up: 63*032b8ab6SVijay Mahadevan ierr = ISLocalToGlobalMappingDestroy(<og);CHKERRQ(ierr); 64*032b8ab6SVijay Mahadevan ierr = PetscFree(gindices);CHKERRQ(ierr); 65*032b8ab6SVijay Mahadevan 66*032b8ab6SVijay Mahadevan } else { 67*032b8ab6SVijay Mahadevan ierr = MatSetBlockSize(*J, dmmoab->bs);CHKERRQ(ierr); 68*032b8ab6SVijay Mahadevan } 69*032b8ab6SVijay Mahadevan 70*032b8ab6SVijay Mahadevan /* set preallocation based on different supported Mat types */ 71*032b8ab6SVijay Mahadevan ierr = MatSeqAIJSetPreallocation(*J, innz, nnz);CHKERRQ(ierr); 72*032b8ab6SVijay Mahadevan ierr = MatMPIAIJSetPreallocation(*J, innz, nnz, ionz, onz);CHKERRQ(ierr); 73*032b8ab6SVijay Mahadevan ierr = MatSeqBAIJSetPreallocation(*J, dmmoab->bs, innz, nnz);CHKERRQ(ierr); 74*032b8ab6SVijay Mahadevan ierr = MatMPIBAIJSetPreallocation(*J, dmmoab->bs, innz, nnz, ionz, onz);CHKERRQ(ierr); 75*032b8ab6SVijay Mahadevan 76*032b8ab6SVijay Mahadevan ierr = MatSetDM(*J, dm);CHKERRQ(ierr); 77*032b8ab6SVijay Mahadevan ierr = PetscFree(nnz);CHKERRQ(ierr); 78*032b8ab6SVijay Mahadevan ierr = PetscFree(onz);CHKERRQ(ierr); 79*032b8ab6SVijay Mahadevan PetscFunctionReturn(0); 80*032b8ab6SVijay Mahadevan } 81*032b8ab6SVijay Mahadevan 82*032b8ab6SVijay Mahadevan 83*032b8ab6SVijay Mahadevan #undef __FUNCT__ 84*032b8ab6SVijay Mahadevan #define __FUNCT__ "DMMoab_Compute_NNZ_From_Connectivity" 85*032b8ab6SVijay Mahadevan PetscErrorCode DMMoab_Compute_NNZ_From_Connectivity(DM dm,PetscInt* innz,PetscInt* nnz,PetscInt* ionz,PetscInt* onz,PetscBool isbaij) 86*032b8ab6SVijay Mahadevan { 87*032b8ab6SVijay Mahadevan PetscErrorCode ierr; 88*032b8ab6SVijay Mahadevan PetscInt i,j,k,nloc,count,vpere,bs,nsize,nghost_found; 89*032b8ab6SVijay Mahadevan DM_Moab *dmmoab = (DM_Moab*)dm->data; 90*032b8ab6SVijay Mahadevan const moab::EntityHandle *connect; 91*032b8ab6SVijay Mahadevan moab::Range *vowned=dmmoab->vowned,*elocal=dmmoab->elocal,*eghost=dmmoab->eghost; 92*032b8ab6SVijay Mahadevan moab::Range::iterator iter; 93*032b8ab6SVijay Mahadevan PetscInt *vertex_ids,firstvtx,dof,offset; 94*032b8ab6SVijay Mahadevan moab::Tag id_tag; 95*032b8ab6SVijay Mahadevan moab::ErrorCode merr; 96*032b8ab6SVijay Mahadevan 97*032b8ab6SVijay Mahadevan PetscFunctionBegin; 98*032b8ab6SVijay Mahadevan bs = dmmoab->bs; 99*032b8ab6SVijay Mahadevan nloc = dmmoab->nloc; 100*032b8ab6SVijay Mahadevan nsize = (isbaij ? nloc:nloc*bs); 101*032b8ab6SVijay Mahadevan 102*032b8ab6SVijay Mahadevan ierr = DMMoabGetLocalToGlobalTag(dm,&id_tag);CHKERRQ(ierr); 103*032b8ab6SVijay Mahadevan merr = dmmoab->mbiface->tag_iterate(id_tag,dmmoab->vlocal->begin(),dmmoab->vlocal->end(), 104*032b8ab6SVijay Mahadevan count,reinterpret_cast<void*&>(vertex_ids));MBERRNM(merr); 105*032b8ab6SVijay Mahadevan 106*032b8ab6SVijay Mahadevan merr = dmmoab->mbiface->tag_get_data(dmmoab->ltog_tag,&(*vowned->begin()),1,&firstvtx);MBERRNM(merr); 107*032b8ab6SVijay Mahadevan 108*032b8ab6SVijay Mahadevan for(iter = elocal->begin(),i=0; iter != elocal->end(); iter++,i++) { 109*032b8ab6SVijay Mahadevan /* Get connectivity information in canonical ordering for the local element */ 110*032b8ab6SVijay Mahadevan merr = dmmoab->mbiface->get_connectivity((*iter),connect,vpere);MBERRNM(merr); 111*032b8ab6SVijay Mahadevan 112*032b8ab6SVijay Mahadevan nghost_found=0; 113*032b8ab6SVijay Mahadevan /* loop over vertices and update the number of connectivity */ 114*032b8ab6SVijay Mahadevan for (j=0;j<vpere;j++) { 115*032b8ab6SVijay Mahadevan moab::Range::const_iterator giter = dmmoab->vghost->find(connect[j]); 116*032b8ab6SVijay Mahadevan if (giter != dmmoab->vghost->end()) nghost_found++; 117*032b8ab6SVijay Mahadevan } 118*032b8ab6SVijay Mahadevan 119*032b8ab6SVijay Mahadevan /* loop over vertices and update the number of connectivity */ 120*032b8ab6SVijay Mahadevan for (j=0;j<vpere;j++) { 121*032b8ab6SVijay Mahadevan merr = dmmoab->mbiface->tag_get_data(dmmoab->ltog_tag,&connect[j],1,&dof);MBERRNM(merr); 122*032b8ab6SVijay Mahadevan 123*032b8ab6SVijay Mahadevan moab::Range::const_iterator giter = dmmoab->vghost->find(connect[j]); 124*032b8ab6SVijay Mahadevan if (giter != dmmoab->vghost->end()) continue; 125*032b8ab6SVijay Mahadevan 126*032b8ab6SVijay Mahadevan /* if block format, then all the block data are local only */ 127*032b8ab6SVijay Mahadevan if (!isbaij) { 128*032b8ab6SVijay Mahadevan offset=(dof-firstvtx)*bs; 129*032b8ab6SVijay Mahadevan for (k=0;k<bs;k++) { 130*032b8ab6SVijay Mahadevan nnz[offset+k]+=vpere-1; 131*032b8ab6SVijay Mahadevan onz[offset+k]+=nghost_found; /* found a ghost non-owned node */ 132*032b8ab6SVijay Mahadevan } 133*032b8ab6SVijay Mahadevan } 134*032b8ab6SVijay Mahadevan else { 135*032b8ab6SVijay Mahadevan nnz[(dof-firstvtx)]+=vpere-1; 136*032b8ab6SVijay Mahadevan onz[(dof-firstvtx)]+=nghost_found; /* found a ghost non-owned node */ 137*032b8ab6SVijay Mahadevan } 138*032b8ab6SVijay Mahadevan } 139*032b8ab6SVijay Mahadevan } 140*032b8ab6SVijay Mahadevan 141*032b8ab6SVijay Mahadevan for(iter = eghost->begin(),i=0; iter != eghost->end(); iter++,i++) { 142*032b8ab6SVijay Mahadevan /* Get connectivity information in canonical ordering for the local element */ 143*032b8ab6SVijay Mahadevan merr = dmmoab->mbiface->get_connectivity((*iter),connect,vpere);MBERRNM(merr); 144*032b8ab6SVijay Mahadevan 145*032b8ab6SVijay Mahadevan nghost_found=0; 146*032b8ab6SVijay Mahadevan /* loop over vertices and update the number of connectivity */ 147*032b8ab6SVijay Mahadevan for (j=0;j<vpere;j++) { 148*032b8ab6SVijay Mahadevan moab::Range::const_iterator giter = dmmoab->vghost->find(connect[j]); 149*032b8ab6SVijay Mahadevan if (giter != dmmoab->vghost->end()) nghost_found++; 150*032b8ab6SVijay Mahadevan } 151*032b8ab6SVijay Mahadevan 152*032b8ab6SVijay Mahadevan /* loop over vertices and update the number of connectivity */ 153*032b8ab6SVijay Mahadevan for (j=0;j<vpere;j++) { 154*032b8ab6SVijay Mahadevan moab::Range::const_iterator giter = dmmoab->vghost->find(connect[j]); 155*032b8ab6SVijay Mahadevan if (giter != dmmoab->vghost->end()) continue; 156*032b8ab6SVijay Mahadevan 157*032b8ab6SVijay Mahadevan merr = dmmoab->mbiface->tag_get_data(dmmoab->ltog_tag,&connect[j],1,&dof);MBERRNM(merr); 158*032b8ab6SVijay 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); 159*032b8ab6SVijay Mahadevan 160*032b8ab6SVijay Mahadevan /* if block format, then all the block data are local only */ 161*032b8ab6SVijay Mahadevan if (!isbaij) { 162*032b8ab6SVijay Mahadevan offset=(dof-firstvtx)*bs; 163*032b8ab6SVijay Mahadevan for (k=0;k<bs;k++) 164*032b8ab6SVijay Mahadevan onz[offset+k]+=nghost_found; 165*032b8ab6SVijay Mahadevan } 166*032b8ab6SVijay Mahadevan else { 167*032b8ab6SVijay Mahadevan onz[(dof-firstvtx)]+=nghost_found; 168*032b8ab6SVijay Mahadevan } 169*032b8ab6SVijay Mahadevan } 170*032b8ab6SVijay Mahadevan } 171*032b8ab6SVijay Mahadevan 172*032b8ab6SVijay Mahadevan if (innz) *innz=0; 173*032b8ab6SVijay Mahadevan if (ionz) *ionz=0; 174*032b8ab6SVijay Mahadevan for (i=0;i<nsize;i++) { 175*032b8ab6SVijay Mahadevan nnz[i]+=1; /* self count the node */ 176*032b8ab6SVijay Mahadevan if (!isbaij) { 177*032b8ab6SVijay Mahadevan nnz[i]*=bs; 178*032b8ab6SVijay Mahadevan onz[i]*=bs; 179*032b8ab6SVijay Mahadevan } 180*032b8ab6SVijay Mahadevan if (innz && nnz[i]>*innz) *innz=nnz[i]; 181*032b8ab6SVijay Mahadevan if (ionz && onz[i]>*ionz) *ionz=onz[i]; 182*032b8ab6SVijay Mahadevan } 183*032b8ab6SVijay Mahadevan // PetscPrintf (PETSC_COMM_WORLD, "\n Maximum NNZ = %D and ONZ = %D \n", *innz,*ionz); 184*032b8ab6SVijay Mahadevan 185*032b8ab6SVijay Mahadevan /* Print for debug purposes only */ 186*032b8ab6SVijay Mahadevan /* 187*032b8ab6SVijay Mahadevan IS nzis,onzis; 188*032b8ab6SVijay Mahadevan ierr = ISCreateGeneral(((PetscObject)dm)->comm, nsize, nnz, PETSC_COPY_VALUES, &nzis);CHKERRQ(ierr); 189*032b8ab6SVijay Mahadevan ierr = ISCreateGeneral(((PetscObject)dm)->comm, nsize, onz, PETSC_COPY_VALUES, &onzis);CHKERRQ(ierr); 190*032b8ab6SVijay Mahadevan ierr = ISView(nzis, 0);CHKERRQ(ierr); 191*032b8ab6SVijay Mahadevan ierr = ISView(onzis, 0);CHKERRQ(ierr); 192*032b8ab6SVijay Mahadevan ierr = ISDestroy(&onzis);CHKERRQ(ierr); 193*032b8ab6SVijay Mahadevan ierr = ISDestroy(&onzis);CHKERRQ(ierr); 194*032b8ab6SVijay Mahadevan std::cin.get(); 195*032b8ab6SVijay Mahadevan */ 196*032b8ab6SVijay Mahadevan 197*032b8ab6SVijay Mahadevan 198*032b8ab6SVijay Mahadevan PetscFunctionReturn(0); 199*032b8ab6SVijay Mahadevan } 200*032b8ab6SVijay Mahadevan 201*032b8ab6SVijay Mahadevan 202