xref: /petsc/src/mat/tests/ex130.c (revision 856bee69f0e0908e75ff837867b1777dfb1ced96)
1 
2 static char help[] = "Tests external direct solvers. Simplified from ex125.c\n\
3 Example: mpiexec -n <np> ./ex130 -f <matrix binary file> -mat_solver_type 1 -mat_superlu_equil \n\n";
4 
5 #include <petscmat.h>
6 
7 int main(int argc, char **args)
8 {
9   Mat           A, F;
10   Vec           u, x, b;
11   PetscMPIInt   rank, size;
12   PetscInt      m, n, nfact, ipack = 0;
13   PetscReal     norm, tol = 1.e-12, Anorm;
14   IS            perm, iperm;
15   MatFactorInfo info;
16   PetscBool     flg, testMatSolve = PETSC_TRUE;
17   PetscViewer   fd;                       /* viewer */
18   char          file[PETSC_MAX_PATH_LEN]; /* input file name */
19 
20   PetscFunctionBeginUser;
21   PetscCall(PetscInitialize(&argc, &args, (char *)0, help));
22   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
23   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
24 
25   /* Determine file from which we read the matrix A */
26   PetscCall(PetscOptionsGetString(NULL, NULL, "-f", file, sizeof(file), &flg));
27   PetscCheck(flg, PETSC_COMM_WORLD, PETSC_ERR_USER, "Must indicate binary file with the -f option");
28 
29   /* Load matrix A */
30   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, file, FILE_MODE_READ, &fd));
31   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
32   PetscCall(MatLoad(A, fd));
33   PetscCall(VecCreate(PETSC_COMM_WORLD, &b));
34   PetscCall(VecLoad(b, fd));
35   PetscCall(PetscViewerDestroy(&fd));
36   PetscCall(MatGetLocalSize(A, &m, &n));
37   PetscCheck(m == n, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "This example is not intended for rectangular matrices (%d, %d)", m, n);
38   PetscCall(MatNorm(A, NORM_INFINITY, &Anorm));
39 
40   /* Create vectors */
41   PetscCall(VecDuplicate(b, &x));
42   PetscCall(VecDuplicate(x, &u)); /* save the true solution */
43 
44   /* Test LU Factorization */
45   PetscCall(MatGetOrdering(A, MATORDERINGNATURAL, &perm, &iperm));
46 
47   PetscCall(PetscOptionsGetInt(NULL, NULL, "-mat_solver_type", &ipack, NULL));
48   switch (ipack) {
49   case 1:
50 #if defined(PETSC_HAVE_SUPERLU)
51     if (rank == 0) printf(" SUPERLU LU:\n");
52     PetscCall(MatGetFactor(A, MATSOLVERSUPERLU, MAT_FACTOR_LU, &F));
53     break;
54 #endif
55   case 2:
56 #if defined(PETSC_HAVE_MUMPS)
57     if (rank == 0) printf(" MUMPS LU:\n");
58     PetscCall(MatGetFactor(A, MATSOLVERMUMPS, MAT_FACTOR_LU, &F));
59     {
60       /* test mumps options */
61       PetscInt icntl_7 = 5;
62       PetscCall(MatMumpsSetIcntl(F, 7, icntl_7));
63     }
64     break;
65 #endif
66   default:
67     if (rank == 0) printf(" PETSC LU:\n");
68     PetscCall(MatGetFactor(A, MATSOLVERPETSC, MAT_FACTOR_LU, &F));
69   }
70 
71   info.fill = 5.0;
72   PetscCall(MatLUFactorSymbolic(F, A, perm, iperm, &info));
73 
74   for (nfact = 0; nfact < 1; nfact++) {
75     if (rank == 0) printf(" %d-the LU numfactorization \n", nfact);
76     PetscCall(MatLUFactorNumeric(F, A, &info));
77 
78     /* Test MatSolve() */
79     if (testMatSolve) {
80       PetscCall(MatSolve(F, b, x));
81 
82       /* Check the residual */
83       PetscCall(MatMult(A, x, u));
84       PetscCall(VecAXPY(u, -1.0, b));
85       PetscCall(VecNorm(u, NORM_INFINITY, &norm));
86       if (norm > tol) {
87         if (rank == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "MatSolve: rel residual %g/%g = %g, LU numfact %d\n", norm, Anorm, norm / Anorm, nfact));
88       }
89     }
90   }
91 
92   /* Free data structures */
93   PetscCall(MatDestroy(&A));
94   PetscCall(MatDestroy(&F));
95   PetscCall(ISDestroy(&perm));
96   PetscCall(ISDestroy(&iperm));
97   PetscCall(VecDestroy(&x));
98   PetscCall(VecDestroy(&b));
99   PetscCall(VecDestroy(&u));
100   PetscCall(PetscFinalize());
101   return 0;
102 }
103