Actual source code: amd.c

petsc-3.13.6 2020-09-29
Report Typos and Errors

  2:  #include <petscmat.h>
  3:  #include <petsc/private/matorderimpl.h>
  4: #include <amd.h>

  6: #if defined(PETSC_USE_64BIT_INDICES)
  7: #  define amd_AMD_defaults amd_l_defaults
  8: /* the type casts are needed because PetscInt is long long while SuiteSparse_long is long and compilers warn even when they are identical */
  9: #  define amd_AMD_order(a,b,c,d,e,f)    amd_l_order((SuiteSparse_long)a,(SuiteSparse_long*)b,(SuiteSparse_long*)c,(SuiteSparse_long*)d,e,f)
 10: #else
 11: #  define amd_AMD_defaults amd_defaults
 12: #  define amd_AMD_order    amd_order
 13: #endif

 15: /*
 16:     MatGetOrdering_AMD - Find the Approximate Minimum Degree ordering

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

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

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

 47:   tval = (PetscBool)Control[AMD_AGGRESSIVE];
 48:   PetscOptionsBool("-mat_ordering_amd_aggressive","use aggressive absorption","None",tval,&tval,NULL);
 49:   Control[AMD_AGGRESSIVE] = (double)tval;

 51:   PetscOptionsEnd();

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

 72:   ISCreateGeneral(PETSC_COMM_SELF,nrow,perm,PETSC_COPY_VALUES,row);
 73:   ISCreateGeneral(PETSC_COMM_SELF,nrow,perm,PETSC_OWN_POINTER,col);
 74:   return(0);
 75: }