Actual source code: bipartite.c
petsc-3.7.3 2016-08-01
1: #include <petsc/private/matimpl.h> /*I "petscmat.h" I*/
2: #include <petscsf.h>
6: PETSC_EXTERN PetscErrorCode MatColoringCreateBipartiteGraph(MatColoring mc,PetscSF *etoc,PetscSF *etor)
7: {
8: PetscErrorCode ierr;
9: PetscInt nentries,ncolentries,idx;
10: PetscInt i,j,rs,re,cs,ce,cn;
11: PetscInt *rowleaf,*colleaf,*rowdata;
12: PetscInt ncol;
13: const PetscScalar *vcol;
14: const PetscInt *icol;
15: const PetscInt *coldegrees,*rowdegrees;
16: Mat m = mc->mat;
19: MatGetOwnershipRange(m,&rs,&re);
20: MatGetOwnershipRangeColumn(m,&cs,&ce);
21: cn = ce-cs;
22: nentries=0;
23: for (i=rs;i<re;i++) {
24: MatGetRow(m,i,&ncol,NULL,&vcol);
25: for (j=0;j<ncol;j++) {
26: nentries++;
27: }
28: MatRestoreRow(m,i,&ncol,NULL,&vcol);
29: }
30: PetscMalloc1(nentries,&rowleaf);
31: PetscMalloc1(nentries,&rowdata);
32: idx=0;
33: for (i=rs;i<re;i++) {
34: MatGetRow(m,i,&ncol,&icol,&vcol);
35: for (j=0;j<ncol;j++) {
36: rowleaf[idx] = icol[j];
37: rowdata[idx] = i;
38: idx++;
39: }
40: MatRestoreRow(m,i,&ncol,&icol,&vcol);
41: }
42: if (idx != nentries) SETERRQ2(PetscObjectComm((PetscObject)m),PETSC_ERR_NOT_CONVERGED,"Bad number of entries %d vs %d",idx,nentries);
43: PetscSFCreate(PetscObjectComm((PetscObject)m),etoc);
44: PetscSFCreate(PetscObjectComm((PetscObject)m),etor);
46: PetscSFSetGraphLayout(*etoc,m->cmap,nentries,NULL,PETSC_COPY_VALUES,rowleaf);
47: PetscSFSetFromOptions(*etoc);
49: /* determine the number of entries in the column matrix */
50: PetscLogEventBegin(Mat_Coloring_Comm,*etoc,0,0,0);
51: PetscSFComputeDegreeBegin(*etoc,&coldegrees);
52: PetscSFComputeDegreeEnd(*etoc,&coldegrees);
53: PetscLogEventEnd(Mat_Coloring_Comm,*etoc,0,0,0);
54: ncolentries=0;
55: for (i=0;i<cn;i++) {
56: ncolentries += coldegrees[i];
57: }
58: PetscMalloc1(ncolentries,&colleaf);
60: /* create the one going the other way by building the leaf set */
61: PetscLogEventBegin(Mat_Coloring_Comm,*etoc,0,0,0);
62: PetscSFGatherBegin(*etoc,MPIU_INT,rowdata,colleaf);
63: PetscSFGatherEnd(*etoc,MPIU_INT,rowdata,colleaf);
64: PetscLogEventEnd(Mat_Coloring_Comm,*etoc,0,0,0);
66: /* this one takes mat entries in *columns* to rows -- you never have to actually be able to order the leaf entries. */
67: PetscSFSetGraphLayout(*etor,m->rmap,ncolentries,NULL,PETSC_COPY_VALUES,colleaf);
68: PetscSFSetFromOptions(*etor);
70: PetscLogEventBegin(Mat_Coloring_Comm,*etor,0,0,0);
71: PetscSFComputeDegreeBegin(*etor,&rowdegrees);
72: PetscSFComputeDegreeEnd(*etor,&rowdegrees);
73: PetscLogEventEnd(Mat_Coloring_Comm,*etor,0,0,0);
75: PetscFree(rowdata);
76: PetscFree(rowleaf);
77: PetscFree(colleaf);
78: return(0);
79: }