Actual source code: amd.c
petsc-3.7.3 2016-08-01
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: */
22: PETSC_INTERN PetscErrorCode MatGetOrdering_AMD(Mat mat,MatOrderingType type,IS *row,IS *col)
23: {
25: PetscInt nrow,*perm;
26: const PetscInt *ia,*ja;
27: int status;
28: PetscReal val;
29: double Control[AMD_CONTROL],Info[AMD_INFO];
30: PetscBool 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_COMM_SELF,PETSC_ERR_SUP,"Cannot get rows for matrix type %s",((PetscObject)mat)->type_name);
40: amd_AMD_defaults(Control);
41: PetscOptionsBegin(PetscObjectComm((PetscObject)mat),((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,NULL);
47: Control[AMD_DENSE] = (double)val;
49: tval = (PetscBool)Control[AMD_AGGRESSIVE];
50: PetscOptionsBool("-mat_ordering_amd_aggressive","use aggressive absorption","None",tval,&tval,NULL);
51: Control[AMD_AGGRESSIVE] = (double)tval;
53: PetscOptionsEnd();
55: PetscMalloc1(nrow,&perm);
56: status = amd_AMD_order(nrow,ia,ja,perm,Control,Info);
57: switch (status) {
58: case AMD_OK: break;
59: case AMD_OK_BUT_JUMBLED:
60: /* The result is fine, but PETSc matrices are supposed to satisfy stricter preconditions, so PETSc considers a
61: * matrix that triggers this error condition to be invalid.
62: */
63: SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_PLIB,"According to AMD, the matrix has unsorted and/or duplicate row indices");
64: case AMD_INVALID:
65: amd_info(Info);
66: SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_PLIB,"According to AMD, the matrix is invalid");
67: case AMD_OUT_OF_MEMORY:
68: SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_MEM,"AMD could not compute ordering");
69: default:
70: SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_LIB,"Unexpected return value");
71: }
72: MatRestoreRowIJ(mat,0,PETSC_FALSE,PETSC_TRUE,NULL,&ia,&ja,&done);
74: ISCreateGeneral(PETSC_COMM_SELF,nrow,perm,PETSC_COPY_VALUES,row);
75: ISCreateGeneral(PETSC_COMM_SELF,nrow,perm,PETSC_OWN_POINTER,col);
76: return(0);
77: }