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 27db66d124SVijay 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 33db66d124SVijay 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) { 48db66d124SVijay 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 60db66d124SVijay 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 */ 69*0df6e276SVijay Mahadevan // ierr = MatSeqAIJSetPreallocation(*J, innz, nnz);CHKERRQ(ierr); 70*0df6e276SVijay Mahadevan // ierr = MatMPIAIJSetPreallocation(*J, innz, nnz, ionz, onz);CHKERRQ(ierr); 71*0df6e276SVijay Mahadevan ierr = MatSeqAIJSetPreallocation(*J, innz, 0);CHKERRQ(ierr); 72*0df6e276SVijay Mahadevan ierr = MatMPIAIJSetPreallocation(*J, innz, 0, ionz, 0);CHKERRQ(ierr); 73032b8ab6SVijay Mahadevan ierr = MatSeqBAIJSetPreallocation(*J, dmmoab->bs, innz, nnz);CHKERRQ(ierr); 74032b8ab6SVijay Mahadevan ierr = MatMPIBAIJSetPreallocation(*J, dmmoab->bs, innz, nnz, ionz, onz);CHKERRQ(ierr); 75032b8ab6SVijay Mahadevan 76032b8ab6SVijay Mahadevan ierr = MatSetDM(*J, dm);CHKERRQ(ierr); 77032b8ab6SVijay Mahadevan ierr = PetscFree(nnz);CHKERRQ(ierr); 78032b8ab6SVijay Mahadevan ierr = PetscFree(onz);CHKERRQ(ierr); 79032b8ab6SVijay Mahadevan PetscFunctionReturn(0); 80032b8ab6SVijay Mahadevan } 81032b8ab6SVijay Mahadevan 82032b8ab6SVijay Mahadevan 83032b8ab6SVijay Mahadevan #undef __FUNCT__ 84032b8ab6SVijay Mahadevan #define __FUNCT__ "DMMoab_Compute_NNZ_From_Connectivity" 85032b8ab6SVijay Mahadevan PetscErrorCode DMMoab_Compute_NNZ_From_Connectivity(DM dm,PetscInt* innz,PetscInt* nnz,PetscInt* ionz,PetscInt* onz,PetscBool isbaij) 86032b8ab6SVijay Mahadevan { 87032b8ab6SVijay Mahadevan PetscErrorCode ierr; 88*0df6e276SVijay Mahadevan PetscInt i,j,f,nloc,vpere,bs,nsize,nghost_found; 89032b8ab6SVijay Mahadevan DM_Moab *dmmoab = (DM_Moab*)dm->data; 90032b8ab6SVijay Mahadevan const moab::EntityHandle *connect; 91*0df6e276SVijay Mahadevan moab::Range *vowned=dmmoab->vowned,*elocal=dmmoab->elocal,*eghost=dmmoab->eghost,adjs,visited; 92032b8ab6SVijay Mahadevan moab::Range::iterator iter; 93*0df6e276SVijay Mahadevan PetscInt dof,doff,ndofs; 94032b8ab6SVijay Mahadevan moab::Tag id_tag; 95032b8ab6SVijay Mahadevan moab::ErrorCode merr; 96*0df6e276SVijay Mahadevan PetscSection section; 97032b8ab6SVijay Mahadevan 98032b8ab6SVijay Mahadevan PetscFunctionBegin; 99032b8ab6SVijay Mahadevan bs = dmmoab->bs; 100032b8ab6SVijay Mahadevan nloc = dmmoab->nloc; 101032b8ab6SVijay Mahadevan nsize = (isbaij ? nloc:nloc*bs); 102032b8ab6SVijay Mahadevan 103032b8ab6SVijay Mahadevan ierr = DMMoabGetLocalToGlobalTag(dm,&id_tag);CHKERRQ(ierr); 104032b8ab6SVijay Mahadevan 105*0df6e276SVijay Mahadevan ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 106032b8ab6SVijay Mahadevan 107032b8ab6SVijay Mahadevan for(iter = elocal->begin(),i=0; iter != elocal->end(); iter++,i++) { 108032b8ab6SVijay Mahadevan /* Get connectivity information in canonical ordering for the local element */ 109032b8ab6SVijay Mahadevan merr = dmmoab->mbiface->get_connectivity((*iter),connect,vpere);MBERRNM(merr); 110032b8ab6SVijay Mahadevan 111032b8ab6SVijay Mahadevan nghost_found=0; 112032b8ab6SVijay Mahadevan /* loop over vertices and update the number of connectivity */ 113032b8ab6SVijay Mahadevan for (j=0;j<vpere;j++) { 114032b8ab6SVijay Mahadevan moab::Range::const_iterator giter = dmmoab->vghost->find(connect[j]); 115032b8ab6SVijay Mahadevan if (giter != dmmoab->vghost->end()) nghost_found++; 116032b8ab6SVijay Mahadevan } 117*0df6e276SVijay Mahadevan // PetscPrintf(PETSC_COMM_WORLD, "Element %D, VPERE %D, nGhost %D\n", i,vpere,nghost_found); 118032b8ab6SVijay Mahadevan 119032b8ab6SVijay Mahadevan /* loop over vertices and update the number of connectivity */ 120032b8ab6SVijay Mahadevan for (j=0;j<vpere;j++) { 121032b8ab6SVijay Mahadevan 122*0df6e276SVijay Mahadevan /* 123*0df6e276SVijay Mahadevan moab::Range::const_iterator giter = visited.find(connect[j]); 124*0df6e276SVijay Mahadevan if (giter != visited.end()) continue; 125032b8ab6SVijay Mahadevan 126*0df6e276SVijay Mahadevan visited.insert(connect[j]); 127*0df6e276SVijay Mahadevan */ 128*0df6e276SVijay Mahadevan 129*0df6e276SVijay Mahadevan merr = dmmoab->mbiface->get_adjacencies(&connect[j],1,1,true,adjs,moab::Interface::INTERSECT); 130*0df6e276SVijay Mahadevan // PetscPrintf(PETSC_COMM_WORLD, "Point %D, nADJS %D\n", connect[j], adjs.size()); 131*0df6e276SVijay Mahadevan 132*0df6e276SVijay Mahadevan merr = dmmoab->mbiface->tag_get_data(id_tag,&connect[j],1,&dof);MBERRNM(merr); 133*0df6e276SVijay Mahadevan ierr = PetscSectionGetDof(section, connect[j], &ndofs); 134*0df6e276SVijay Mahadevan 135*0df6e276SVijay Mahadevan for (f=0;f<dmmoab->nfields;f++) { 136*0df6e276SVijay Mahadevan ierr = PetscSectionGetFieldOffset(section, connect[j], f, &doff); 137*0df6e276SVijay Mahadevan 138*0df6e276SVijay Mahadevan // nnz[doff]+=(vpere-nghost_found-1); /* leave out self to avoid repeats -> node shared by multiple elements */ 139*0df6e276SVijay Mahadevan nnz[doff]+=adjs.size(); /* leave out self to avoid repeats -> node shared by multiple elements */ 140*0df6e276SVijay Mahadevan onz[doff]+=nghost_found; /* add ghost non-owned nodes */ 141*0df6e276SVijay Mahadevan 142*0df6e276SVijay Mahadevan // PetscPrintf(PETSC_COMM_WORLD, "Point %D, Field %D, DOFPP %D \t DOF %D, OFFSET %D : NNZ = %D, ONZ = %D\n", connect[j], f, ndofs, dof, doff, nnz[doff], onz[doff]); 143*0df6e276SVijay Mahadevan adjs.clear(); 144*0df6e276SVijay Mahadevan } 145*0df6e276SVijay Mahadevan } 146*0df6e276SVijay Mahadevan } 147*0df6e276SVijay Mahadevan 148*0df6e276SVijay Mahadevan visited.clear(); 149*0df6e276SVijay Mahadevan 150*0df6e276SVijay Mahadevan if (innz) *innz=0; 151*0df6e276SVijay Mahadevan if (ionz) *ionz=0; 152*0df6e276SVijay Mahadevan for (i=0;i<nsize;i++) { 153*0df6e276SVijay Mahadevan nnz[i]+=1; /* self count the node */ 154032b8ab6SVijay Mahadevan if (!isbaij) { 155*0df6e276SVijay Mahadevan nnz[i]*=bs; 156*0df6e276SVijay Mahadevan onz[i]*=bs; 157032b8ab6SVijay Mahadevan } 158*0df6e276SVijay Mahadevan // PetscPrintf(PETSC_COMM_WORLD, "Index %D: NNZ = %D, ONZ = %D\n", i, nnz[i], onz[i]); 159*0df6e276SVijay Mahadevan 160*0df6e276SVijay Mahadevan if (innz && nnz[i]>*innz) *innz=nnz[i]; 161*0df6e276SVijay Mahadevan if (ionz && onz[i]>*ionz) *ionz=onz[i]; 162032b8ab6SVijay Mahadevan } 163*0df6e276SVijay Mahadevan // PetscPrintf(PETSC_COMM_WORLD, "MAX: NNZ = %D, ONZ = %D\n", *innz, *ionz); 164*0df6e276SVijay Mahadevan PetscFunctionReturn(0); 165032b8ab6SVijay Mahadevan } 166032b8ab6SVijay Mahadevan 167*0df6e276SVijay Mahadevan 168*0df6e276SVijay Mahadevan #if 0 169*0df6e276SVijay Mahadevan 170*0df6e276SVijay Mahadevan #undef __FUNCT__ 171*0df6e276SVijay Mahadevan #define __FUNCT__ "DMMoab_Compute_NNZ_From_Connectivity_OLD" 172*0df6e276SVijay Mahadevan PetscErrorCode DMMoab_Compute_NNZ_From_Connectivity_OLD(DM dm,PetscInt* innz,PetscInt* nnz,PetscInt* ionz,PetscInt* onz,PetscBool isbaij) 173*0df6e276SVijay Mahadevan { 174*0df6e276SVijay Mahadevan PetscErrorCode ierr; 175*0df6e276SVijay Mahadevan PetscInt i,j,k,f,nloc,count,vpere,bs,nsize,nghost_found; 176*0df6e276SVijay Mahadevan DM_Moab *dmmoab = (DM_Moab*)dm->data; 177*0df6e276SVijay Mahadevan const moab::EntityHandle *connect; 178*0df6e276SVijay Mahadevan moab::Range *vowned=dmmoab->vowned,*elocal=dmmoab->elocal,*eghost=dmmoab->eghost,adjs; 179*0df6e276SVijay Mahadevan moab::Range::iterator iter; 180*0df6e276SVijay Mahadevan PetscInt *vertex_ids,firstvtx,dof,doff,ndofs,offset; 181*0df6e276SVijay Mahadevan moab::Tag id_tag; 182*0df6e276SVijay Mahadevan std::vector<moab::EntityHandle> storage; 183*0df6e276SVijay Mahadevan moab::ErrorCode merr; 184*0df6e276SVijay Mahadevan PetscSection section; 185*0df6e276SVijay Mahadevan 186*0df6e276SVijay Mahadevan PetscFunctionBegin; 187*0df6e276SVijay Mahadevan bs = dmmoab->bs; 188*0df6e276SVijay Mahadevan nloc = dmmoab->nloc; 189*0df6e276SVijay Mahadevan nsize = (isbaij ? nloc:nloc*bs); 190*0df6e276SVijay Mahadevan 191*0df6e276SVijay Mahadevan ierr = DMMoabGetLocalToGlobalTag(dm,&id_tag);CHKERRQ(ierr); 192*0df6e276SVijay Mahadevan 193*0df6e276SVijay Mahadevan merr = dmmoab->mbiface->tag_get_data(id_tag,&(*vowned->begin()),1,&firstvtx);MBERRNM(merr); 194*0df6e276SVijay Mahadevan //firstvtx = vertex_ids[0]; 195*0df6e276SVijay Mahadevan 196*0df6e276SVijay Mahadevan ierr = DMGetDefaultSection(dm, §ion);CHKERRQ(ierr); 197*0df6e276SVijay Mahadevan 198*0df6e276SVijay Mahadevan for(iter = elocal->begin(),i=0; iter != elocal->end(); iter++,i++) { 199032b8ab6SVijay Mahadevan /* Get connectivity information in canonical ordering for the local element */ 200*0df6e276SVijay Mahadevan merr = dmmoab->mbiface->get_connectivity((*iter),connect,vpere,false,&storage);MBERRNM(merr); 201032b8ab6SVijay Mahadevan 202032b8ab6SVijay Mahadevan nghost_found=0; 203032b8ab6SVijay Mahadevan /* loop over vertices and update the number of connectivity */ 204032b8ab6SVijay Mahadevan for (j=0;j<vpere;j++) { 205032b8ab6SVijay Mahadevan moab::Range::const_iterator giter = dmmoab->vghost->find(connect[j]); 206032b8ab6SVijay Mahadevan if (giter != dmmoab->vghost->end()) nghost_found++; 207032b8ab6SVijay Mahadevan } 208032b8ab6SVijay Mahadevan 209032b8ab6SVijay Mahadevan /* loop over vertices and update the number of connectivity */ 210032b8ab6SVijay Mahadevan for (j=0;j<vpere;j++) { 211*0df6e276SVijay Mahadevan 212*0df6e276SVijay Mahadevan merr = dmmoab->mbiface->get_adjacencies(&connect[j],1,1,true,adjs,moab::Interface::INTERSECT); 213*0df6e276SVijay Mahadevan PetscPrintf(PETSC_COMM_WORLD, "Point %D, nADJS %D\n", connect[j], adjs.size()); 214*0df6e276SVijay Mahadevan 215032b8ab6SVijay Mahadevan moab::Range::const_iterator giter = dmmoab->vghost->find(connect[j]); 216032b8ab6SVijay Mahadevan if (giter != dmmoab->vghost->end()) continue; 217032b8ab6SVijay Mahadevan 218*0df6e276SVijay Mahadevan merr = dmmoab->mbiface->tag_get_data(id_tag,&connect[j],1,&dof);MBERRNM(merr); 219*0df6e276SVijay Mahadevan ierr = PetscSectionGetDof(section, connect[j], &ndofs); 220032b8ab6SVijay Mahadevan 221*0df6e276SVijay Mahadevan for (f=0;f<dmmoab->nfields;f++) { 222*0df6e276SVijay Mahadevan ierr = PetscSectionGetFieldOffset(section, connect[j], f, &doff); 223*0df6e276SVijay Mahadevan 224*0df6e276SVijay Mahadevan // nnz[doff]+=ndofs*(vpere-nghost_found-1); /* leave out self to avoid repeats -> node shared by multiple elements */ 225*0df6e276SVijay Mahadevan nnz[doff]+=adjs.size()-nghost_found-1; /* leave out self to avoid repeats -> node shared by multiple elements */ 226*0df6e276SVijay Mahadevan onz[doff]+=nghost_found; /* add ghost non-owned nodes */ 227*0df6e276SVijay Mahadevan 228*0df6e276SVijay Mahadevan PetscPrintf(PETSC_COMM_WORLD, "Point %D, Field %D, DOFPP %D \t DOF %D, OFFSET %D : NNZ = %D, ONZ = %D\n", connect[j], f, ndofs, dof, doff, nnz[doff], onz[doff]); 229*0df6e276SVijay Mahadevan adjs.clear(); 230032b8ab6SVijay Mahadevan } 231*0df6e276SVijay Mahadevan 232*0df6e276SVijay Mahadevan // /* if block format, then all the block data are local only */ 233*0df6e276SVijay Mahadevan // if (!isbaij) { 234*0df6e276SVijay Mahadevan // offset=(dof-firstvtx)*bs; 235*0df6e276SVijay Mahadevan // for (k=0;k<bs;k++) { 236*0df6e276SVijay Mahadevan // nnz[offset+k]+=vpere-nghost_found-1; /* leave out self to avoid repeats -> node shared by multiple elements */ 237*0df6e276SVijay Mahadevan // onz[offset+k]+=nghost_found; /* found a ghost non-owned node */ 238*0df6e276SVijay Mahadevan // } 239*0df6e276SVijay Mahadevan // } 240*0df6e276SVijay Mahadevan // else { 241*0df6e276SVijay Mahadevan // nnz[(dof-firstvtx)]+=vpere-nghost_found-1; /* leave out self to avoid repeats -> node shared by multiple elements */ 242*0df6e276SVijay Mahadevan // onz[(dof-firstvtx)]+=nghost_found; /* found a ghost non-owned node */ 243*0df6e276SVijay Mahadevan // } 244*0df6e276SVijay Mahadevan 245032b8ab6SVijay Mahadevan } 246032b8ab6SVijay Mahadevan } 247*0df6e276SVijay Mahadevan 248*0df6e276SVijay Mahadevan /* loop through ghosted elements that contain non-ghosted (locally owned) vertices */ 249*0df6e276SVijay Mahadevan for(iter = eghost->begin(),i=0; iter != eghost->end(); iter++,i++) { 250*0df6e276SVijay Mahadevan /* Get connectivity information in canonical ordering for the local element */ 251*0df6e276SVijay Mahadevan merr = dmmoab->mbiface->get_connectivity((*iter),connect,vpere,false,&storage);MBERRNM(merr); 252*0df6e276SVijay Mahadevan 253*0df6e276SVijay Mahadevan nghost_found=0; 254*0df6e276SVijay Mahadevan /* loop over vertices and update the number of connectivity */ 255*0df6e276SVijay Mahadevan for (j=0;j<vpere;j++) { 256*0df6e276SVijay Mahadevan moab::Range::const_iterator giter = dmmoab->vghost->find(connect[j]); 257*0df6e276SVijay Mahadevan if (giter != dmmoab->vghost->end()) nghost_found++; 258*0df6e276SVijay Mahadevan } 259*0df6e276SVijay Mahadevan 260*0df6e276SVijay Mahadevan if (nghost_found == vpere) continue; /* all vertices are ghosted.. */ 261*0df6e276SVijay Mahadevan 262*0df6e276SVijay Mahadevan /* loop over vertices and update the number of connectivity */ 263*0df6e276SVijay Mahadevan for (j=0;j<vpere;j++) { 264*0df6e276SVijay Mahadevan moab::Range::const_iterator giter = dmmoab->vghost->find(connect[j]); 265*0df6e276SVijay Mahadevan if (giter != dmmoab->vghost->end()) continue; 266*0df6e276SVijay Mahadevan 267*0df6e276SVijay Mahadevan merr = dmmoab->mbiface->tag_get_data(id_tag,&connect[j],1,&dof);MBERRNM(merr); 268*0df6e276SVijay Mahadevan ierr = PetscSectionGetDof(section, connect[j], &ndofs); 269*0df6e276SVijay Mahadevan 270*0df6e276SVijay Mahadevan for (f=0;f<dmmoab->nfields;f++) { 271*0df6e276SVijay Mahadevan ierr = PetscSectionGetFieldOffset(section, connect[j], f, &doff); 272*0df6e276SVijay Mahadevan 273*0df6e276SVijay Mahadevan onz[doff]+=nghost_found; /* add ghost non-owned nodes */ 274*0df6e276SVijay Mahadevan 275*0df6e276SVijay Mahadevan PetscPrintf(PETSC_COMM_WORLD, "Point %D, Field %D, DOF %D, OFFSET %D : ONZ = %D\n", connect[j], f, dof, doff, onz[doff]); 276*0df6e276SVijay Mahadevan 277*0df6e276SVijay Mahadevan } 278*0df6e276SVijay Mahadevan 279*0df6e276SVijay 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); 280*0df6e276SVijay Mahadevan // 281*0df6e276SVijay Mahadevan // /* if block format, then all the block data are local only */ 282*0df6e276SVijay Mahadevan // if (!isbaij) { 283*0df6e276SVijay Mahadevan // offset=(dof-firstvtx)*bs; 284*0df6e276SVijay Mahadevan // for (k=0;k<bs;k++) 285*0df6e276SVijay Mahadevan // onz[offset+k]+=nghost_found; 286*0df6e276SVijay Mahadevan // } 287*0df6e276SVijay Mahadevan // else { 288*0df6e276SVijay Mahadevan // onz[(dof-firstvtx)]+=nghost_found; 289*0df6e276SVijay Mahadevan // } 290*0df6e276SVijay Mahadevan } 291032b8ab6SVijay Mahadevan } 292032b8ab6SVijay Mahadevan 293032b8ab6SVijay Mahadevan if (innz) *innz=0; 294032b8ab6SVijay Mahadevan if (ionz) *ionz=0; 295032b8ab6SVijay Mahadevan for (i=0;i<nsize;i++) { 296032b8ab6SVijay Mahadevan nnz[i]+=1; /* self count the node */ 297032b8ab6SVijay Mahadevan if (!isbaij) { 298032b8ab6SVijay Mahadevan nnz[i]*=bs; 299032b8ab6SVijay Mahadevan onz[i]*=bs; 300032b8ab6SVijay Mahadevan } 301*0df6e276SVijay Mahadevan PetscPrintf(PETSC_COMM_WORLD, "Index %D: NNZ = %D, ONZ = %D\n", i, nnz[i], onz[i]); 302*0df6e276SVijay Mahadevan 303032b8ab6SVijay Mahadevan if (innz && nnz[i]>*innz) *innz=nnz[i]; 304032b8ab6SVijay Mahadevan if (ionz && onz[i]>*ionz) *ionz=onz[i]; 305032b8ab6SVijay Mahadevan } 306*0df6e276SVijay Mahadevan PetscPrintf(PETSC_COMM_WORLD, "MAX: NNZ = %D, ONZ = %D\n", *innz, *ionz); 307032b8ab6SVijay Mahadevan PetscFunctionReturn(0); 308032b8ab6SVijay Mahadevan } 309032b8ab6SVijay Mahadevan 310*0df6e276SVijay Mahadevan #endif 311032b8ab6SVijay Mahadevan 312