Actual source code: bipartite.c

petsc-3.8.4 2018-03-24
Report Typos and Errors
  1:  #include <petsc/private/matimpl.h>
  2:  #include <petscsf.h>

  4: PETSC_EXTERN PetscErrorCode MatColoringCreateBipartiteGraph(MatColoring mc,PetscSF *etoc,PetscSF *etor)
  5: {
  6:   PetscErrorCode    ierr;
  7:   PetscInt          nentries,ncolentries,idx;
  8:   PetscInt          i,j,rs,re,cs,ce,cn;
  9:   PetscInt          *rowleaf,*colleaf,*rowdata;
 10:   PetscInt          ncol;
 11:   const PetscScalar *vcol;
 12:   const PetscInt    *icol;
 13:   const PetscInt    *coldegrees,*rowdegrees;
 14:   Mat               m = mc->mat;

 17:   MatGetOwnershipRange(m,&rs,&re);
 18:   MatGetOwnershipRangeColumn(m,&cs,&ce);
 19:   cn = ce-cs;
 20:   nentries=0;
 21:   for (i=rs;i<re;i++) {
 22:     MatGetRow(m,i,&ncol,NULL,&vcol);
 23:     for (j=0;j<ncol;j++) {
 24:       nentries++;
 25:     }
 26:     MatRestoreRow(m,i,&ncol,NULL,&vcol);
 27:   }
 28:   PetscMalloc1(nentries,&rowleaf);
 29:   PetscMalloc1(nentries,&rowdata);
 30:   idx=0;
 31:   for (i=rs;i<re;i++) {
 32:     MatGetRow(m,i,&ncol,&icol,&vcol);
 33:     for (j=0;j<ncol;j++) {
 34:       rowleaf[idx] = icol[j];
 35:       rowdata[idx] = i;
 36:       idx++;
 37:     }
 38:     MatRestoreRow(m,i,&ncol,&icol,&vcol);
 39:   }
 40:   if (idx != nentries) SETERRQ2(PetscObjectComm((PetscObject)m),PETSC_ERR_NOT_CONVERGED,"Bad number of entries %d vs %d",idx,nentries);
 41:   PetscSFCreate(PetscObjectComm((PetscObject)m),etoc);
 42:   PetscSFCreate(PetscObjectComm((PetscObject)m),etor);

 44:   PetscSFSetGraphLayout(*etoc,m->cmap,nentries,NULL,PETSC_COPY_VALUES,rowleaf);
 45:   PetscSFSetFromOptions(*etoc);

 47:   /* determine the number of entries in the column matrix */
 48:   PetscLogEventBegin(MATCOLORING_Comm,*etoc,0,0,0);
 49:   PetscSFComputeDegreeBegin(*etoc,&coldegrees);
 50:   PetscSFComputeDegreeEnd(*etoc,&coldegrees);
 51:   PetscLogEventEnd(MATCOLORING_Comm,*etoc,0,0,0);
 52:   ncolentries=0;
 53:   for (i=0;i<cn;i++) {
 54:     ncolentries += coldegrees[i];
 55:   }
 56:   PetscMalloc1(ncolentries,&colleaf);

 58:   /* create the one going the other way by building the leaf set */
 59:   PetscLogEventBegin(MATCOLORING_Comm,*etoc,0,0,0);
 60:   PetscSFGatherBegin(*etoc,MPIU_INT,rowdata,colleaf);
 61:   PetscSFGatherEnd(*etoc,MPIU_INT,rowdata,colleaf);
 62:   PetscLogEventEnd(MATCOLORING_Comm,*etoc,0,0,0);

 64:   /* this one takes mat entries in *columns* to rows -- you never have to actually be able to order the leaf entries. */
 65:   PetscSFSetGraphLayout(*etor,m->rmap,ncolentries,NULL,PETSC_COPY_VALUES,colleaf);
 66:   PetscSFSetFromOptions(*etor);

 68:   PetscLogEventBegin(MATCOLORING_Comm,*etor,0,0,0);
 69:   PetscSFComputeDegreeBegin(*etor,&rowdegrees);
 70:   PetscSFComputeDegreeEnd(*etor,&rowdegrees);
 71:   PetscLogEventEnd(MATCOLORING_Comm,*etor,0,0,0);

 73:   PetscFree(rowdata);
 74:   PetscFree(rowleaf);
 75:   PetscFree(colleaf);
 76:   return(0);
 77: }