Actual source code: amd.c

  1: #define PETSCMAT_DLL

 3:  #include petscmat.h
 4:  #include ../src/mat/order/order.h
  5: #include <amd.h>

  7: #if defined(PETSC_USE_64BIT_INDICES)
  8: #  define amd_AMD_defaults amd_l_defaults
  9: #  define amd_AMD_order    amd_l_order
 10: #else
 11: #  define amd_AMD_defaults amd_defaults
 12: #  define amd_AMD_order    amd_order
 13: #endif

 16: /*
 17:     MatOrdering_AMD - Find the Approximate Minimum Degree ordering

 19:     This provides an interface to Tim Davis' AMD package (used by UMFPACK, CHOLMOD, Matlab, etc).
 20: */
 23: PetscErrorCode  MatOrdering_AMD(Mat mat,const MatOrderingType type,IS *row,IS *col)
 24: {
 26:   PetscInt       nrow,*ia,*ja,*perm;
 27:   int            status;
 28:   PetscReal      val;
 29:   double         Control[AMD_CONTROL],Info[AMD_INFO];
 30:   PetscTruth     tval,done;

 33:   /*
 34:      AMD does not require that the matrix be symmetric (it does so internally,
 35:      at least in so far as computing orderings for A+A^T.
 36:   */
 37:   MatGetRowIJ(mat,0,PETSC_FALSE,PETSC_TRUE,&nrow,&ia,&ja,&done);
 38:   if (!done) SETERRQ1(PETSC_ERR_SUP,"Cannot get rows for matrix type %s",((PetscObject)mat)->type_name);

 40:   amd_AMD_defaults(Control);
 41:   PetscOptionsBegin(((PetscObject)mat)->comm,((PetscObject)mat)->prefix,"AMD Options","Mat");
 42:   /*
 43:     We have to use temporary values here because AMD always uses double, even though PetscReal may be single
 44:   */
 45:   val = (PetscReal)Control[AMD_DENSE];
 46:   PetscOptionsReal("-mat_ordering_amd_dense","threshold for \"dense\" rows/columns","None",val,&val,PETSC_NULL);
 47:   Control[AMD_DENSE] = (double)val;

 49:   tval = (PetscTruth)Control[AMD_AGGRESSIVE];
 50:   PetscOptionsTruth("-mat_ordering_amd_aggressive","use aggressive absorption","None",tval,&tval,PETSC_NULL);
 51:   Control[AMD_AGGRESSIVE] = (double)tval;
 52:   PetscOptionsEnd();

 54:   PetscMalloc(nrow*sizeof(PetscInt),&perm);
 55:   status = amd_AMD_order(nrow,ia,ja,perm,Control,Info);
 56:   switch (status) {
 57:     case AMD_OK: break;
 58:     case AMD_OK_BUT_JUMBLED:
 59:       /* The result is fine, but PETSc matrices are supposed to satisfy stricter preconditions, so PETSc considers a
 60:       * matrix that triggers this error condition to be invalid.
 61:       */
 62:       SETERRQ(PETSC_ERR_PLIB,"According to AMD, the matrix has unsorted and/or duplicate row indices");
 63:     case AMD_INVALID:
 64:       amd_info(Info);
 65:       SETERRQ(PETSC_ERR_PLIB,"According to AMD, the matrix is invalid");
 66:     case AMD_OUT_OF_MEMORY:
 67:       SETERRQ(PETSC_ERR_MEM,"AMD could not compute ordering");
 68:     default:
 69:       SETERRQ(PETSC_ERR_LIB,"Unexpected return value");
 70:   }
 71:   MatRestoreRowIJ(mat,0,PETSC_FALSE,PETSC_TRUE,&nrow,&ia,&ja,&done);

 73:   ISCreateGeneral(PETSC_COMM_SELF,nrow,perm,row);
 74:   ISCreateGeneral(PETSC_COMM_SELF,nrow,perm,col);
 75:   PetscFree(perm);
 76:   return(0);
 77: }