1 #include <petscmat.h>
2 #include <petsc/private/matorderimpl.h>
3 #include <amd.h>
4
5 #if defined(PETSC_USE_64BIT_INDICES)
6 #define amd_AMD_defaults amd_l_defaults
7 /* the type casts are needed because PetscInt is long long while SuiteSparse_long is long and compilers warn even when they are identical */
8 #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)
9 #else
10 #define amd_AMD_defaults amd_defaults
11 #define amd_AMD_order amd_order
12 #endif
13
14 /*
15 MatGetOrdering_AMD - Find the Approximate Minimum Degree ordering
16
17 This provides an interface to Tim Davis' AMD package (used by UMFPACK, CHOLMOD, MATLAB, etc).
18 */
MatGetOrdering_AMD(Mat mat,MatOrderingType type,IS * row,IS * col)19 PETSC_INTERN PetscErrorCode MatGetOrdering_AMD(Mat mat, MatOrderingType type, IS *row, IS *col)
20 {
21 PetscInt nrow, *perm;
22 const PetscInt *ia, *ja;
23 int status;
24 PetscReal val;
25 double Control[AMD_CONTROL], Info[AMD_INFO];
26 PetscBool tval, done;
27
28 PetscFunctionBegin;
29 /*
30 AMD does not require that the matrix be symmetric (it does so internally,
31 at least in so far as computing orderings for A+A^T.
32 */
33 PetscCall(MatGetRowIJ(mat, 0, PETSC_FALSE, PETSC_TRUE, &nrow, &ia, &ja, &done));
34 PetscCheck(done, PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot get rows for matrix type %s", ((PetscObject)mat)->type_name);
35
36 amd_AMD_defaults(Control);
37 PetscOptionsBegin(PetscObjectComm((PetscObject)mat), ((PetscObject)mat)->prefix, "AMD Options", "Mat");
38 /*
39 We have to use temporary values here because AMD always uses double, even though PetscReal may be single
40 */
41 val = (PetscReal)Control[AMD_DENSE];
42 PetscCall(PetscOptionsReal("-mat_ordering_amd_dense", "threshold for \"dense\" rows/columns", "None", val, &val, NULL));
43 Control[AMD_DENSE] = (double)val;
44
45 tval = (PetscBool)Control[AMD_AGGRESSIVE];
46 PetscCall(PetscOptionsBool("-mat_ordering_amd_aggressive", "use aggressive absorption", "None", tval, &tval, NULL));
47 Control[AMD_AGGRESSIVE] = (double)tval;
48
49 PetscOptionsEnd();
50
51 PetscCall(PetscMalloc1(nrow, &perm));
52 status = amd_AMD_order(nrow, ia, ja, perm, Control, Info);
53 switch (status) {
54 case AMD_OK:
55 break;
56 case AMD_OK_BUT_JUMBLED:
57 /* The result is fine, but PETSc matrices are supposed to satisfy stricter preconditions, so PETSc considers a
58 * matrix that triggers this error condition to be invalid.
59 */
60 SETERRQ(PetscObjectComm((PetscObject)mat), PETSC_ERR_PLIB, "According to AMD, the matrix has unsorted and/or duplicate row indices");
61 case AMD_INVALID:
62 amd_info(Info);
63 SETERRQ(PetscObjectComm((PetscObject)mat), PETSC_ERR_PLIB, "According to AMD, the matrix is invalid");
64 case AMD_OUT_OF_MEMORY:
65 SETERRQ(PetscObjectComm((PetscObject)mat), PETSC_ERR_MEM, "AMD could not compute ordering");
66 default:
67 SETERRQ(PetscObjectComm((PetscObject)mat), PETSC_ERR_LIB, "Unexpected return value");
68 }
69 PetscCall(MatRestoreRowIJ(mat, 0, PETSC_FALSE, PETSC_TRUE, NULL, &ia, &ja, &done));
70
71 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nrow, perm, PETSC_COPY_VALUES, row));
72 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nrow, perm, PETSC_OWN_POINTER, col));
73 PetscFunctionReturn(PETSC_SUCCESS);
74 }
75