Actual source code: wbm.c
petsc-3.5.4 2015-05-23
1: #include <petscmat.h>
2: #include <petsc-private/matorderimpl.h>
4: #if defined(PETSC_HAVE_SUPERLU_DIST)
6: /* SuperLU_DIST bundles f2ced mc64ad_() from HSL */
8: /*
9: SuperLU_dist uses a common flag for both Fortran mangling and BLAS/LAPACK mangling which
10: corresponds to the PETSc BLAS/LAPACK mangling flag (we pass this flag to configure SuperLU_dist)
11: */
12: # if defined(PETSC_BLASLAPACK_CAPS)
13: # define mc64id_ MC64ID
14: # define mc64ad_ MC64AD
15: # elif !defined(PETSC_BLASLAPACK_UNDERSCORE)
16: # define mc64id_ mc64id
17: # define mc64ad_ mc64ad
18: # endif
20: PETSC_EXTERN PetscInt mc64id_(PetscInt*);
21: PETSC_EXTERN PetscInt mc64ad_(const PetscInt*, PetscInt*, PetscInt*, const PetscInt*, const PetscInt*n, PetscScalar*, PetscInt*,
22: PetscInt*, PetscInt*, PetscInt*, PetscInt*, PetscScalar*, PetscInt*, PetscInt*);
23: #endif
25: /*
26: MatGetOrdering_WBM - Find the nonsymmetric reordering of the graph which maximizes the product of diagonal entries,
27: using weighted bipartite graph matching. This is MC64 in the Harwell-Boeing library.
28: */
31: PETSC_EXTERN PetscErrorCode MatGetOrdering_WBM(Mat mat, MatOrderingType type, IS *row, IS *col)
32: {
33: PetscScalar *a, *dw;
34: const PetscInt *ia, *ja;
35: const PetscInt job = 5;
36: PetscInt *perm, nrow, ncol, nnz, liw, *iw, ldw, i;
37: PetscBool done;
38: PetscErrorCode ierr;
41: MatGetRowIJ(mat,1,PETSC_TRUE,PETSC_TRUE,&nrow,&ia,&ja,&done);
42: ncol = nrow;
43: nnz = ia[nrow];
44: if (!done) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"Cannot get rows for matrix");
45: MatSeqAIJGetArray(mat, &a);
46: switch (job) {
47: case 1: liw = 4*nrow + ncol; ldw = 0;break;
48: case 2: liw = 2*nrow + 2*ncol; ldw = ncol;break;
49: case 3: liw = 8*nrow + 2*ncol + nnz; ldw = nnz;break;
50: case 4: liw = 3*nrow + 2*ncol; ldw = 2*ncol + nnz;break;
51: case 5: liw = 3*nrow + 2*ncol; ldw = nrow + 2*ncol + nnz;break;
52: }
54: PetscMalloc3(liw,&iw,ldw,&dw,nrow,&perm);
55: #if defined(PETSC_HAVE_SUPERLU_DIST)
56: {
57: PetscInt num, info[10], icntl[10];
59: mc64id_(icntl);
60: if (ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"HSL mc64id_ returned %d\n",ierr);
61: icntl[0] = 0; /* allow printing error messages (f2c'd code uses if non-negative, ignores value otherwise) */
62: icntl[1] = -1; /* suppress warnings */
63: icntl[2] = -1; /* ignore diagnostic output [default] */
64: icntl[3] = 0; /* perform consistency checks [default] */
65: mc64ad_(&job, &nrow, &nnz, ia, ja, a, &num, perm, &liw, iw, &ldw, dw, icntl, info);
66: if (ierr) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_LIB,"HSL mc64ad_ returned %d\n",ierr);
67: }
68: #else
69: SETERRQ(PetscObjectComm((PetscObject) mat), PETSC_ERR_SUP, "WBM using MC64 does not support complex numbers");
70: #endif
71: MatRestoreRowIJ(mat, 1, PETSC_TRUE, PETSC_TRUE, NULL, &ia, &ja, &done);
72: for (i = 0; i < nrow; ++i) perm[i]--;
73: /* If job == 5, dw[0..ncols] contains the column scaling and dw[ncols..ncols+nrows] contains the row scaling */
74: ISCreateStride(PETSC_COMM_SELF, nrow, 0, 1, row);
75: ISCreateGeneral(PETSC_COMM_SELF,nrow,perm,PETSC_COPY_VALUES,col);
76: PetscFree3(iw,dw,perm);
77: return(0);
78: }