Actual source code: amd.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  2: #include <petscmat.h>
  3: #include <petsc-private/matorderimpl.h>
  4: #define UF_long long long
  5: #define UF_long_max LONG_LONG_MAX
  6: #define UF_long_id "%lld"
  7: #include <amd.h>

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

 17: /*
 18:     MatGetOrdering_AMD - Find the Approximate Minimum Degree ordering

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

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

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

 50:   Control[AMD_DENSE] = (double)val;

 52:   tval = (PetscBool)Control[AMD_AGGRESSIVE];
 53:   PetscOptionsBool("-mat_ordering_amd_aggressive","use aggressive absorption","None",tval,&tval,NULL);

 55:   Control[AMD_AGGRESSIVE] = (double)tval;

 57:   PetscOptionsEnd();

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

 78:   ISCreateGeneral(PETSC_COMM_SELF,nrow,perm,PETSC_COPY_VALUES,row);
 79:   ISCreateGeneral(PETSC_COMM_SELF,nrow,perm,PETSC_OWN_POINTER,col);
 80:   return(0);
 81: }