xref: /petsc/src/mat/tests/ex130.c (revision 39a651e27bfa10b9060071f09b59ca0e2afa5359)
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   PetscErrorCode ierr;
12   PetscMPIInt    rank,size;
13   PetscInt       m,n,nfact,ipack=0;
14   PetscReal      norm,tol=1.e-12,Anorm;
15   IS             perm,iperm;
16   MatFactorInfo  info;
17   PetscBool      flg,testMatSolve=PETSC_TRUE;
18   PetscViewer    fd;              /* viewer */
19   char           file[PETSC_MAX_PATH_LEN]; /* input file name */
20 
21   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
22   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
23   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
24 
25   /* Determine file from which we read the matrix A */
26   CHKERRQ(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   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd));
31   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
32   CHKERRQ(MatLoad(A,fd));
33   CHKERRQ(VecCreate(PETSC_COMM_WORLD,&b));
34   CHKERRQ(VecLoad(b,fd));
35   CHKERRQ(PetscViewerDestroy(&fd));
36   CHKERRQ(MatGetLocalSize(A,&m,&n));
37   PetscCheckFalse(m != n,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ, "This example is not intended for rectangular matrices (%d, %d)", m, n);
38   CHKERRQ(MatNorm(A,NORM_INFINITY,&Anorm));
39 
40   /* Create vectors */
41   CHKERRQ(VecDuplicate(b,&x));
42   CHKERRQ(VecDuplicate(x,&u)); /* save the true solution */
43 
44   /* Test LU Factorization */
45   CHKERRQ(MatGetOrdering(A,MATORDERINGNATURAL,&perm,&iperm));
46 
47   CHKERRQ(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     CHKERRQ(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     CHKERRQ(MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_LU,&F));
59     {
60       /* test mumps options */
61       PetscInt icntl_7 = 5;
62       CHKERRQ(MatMumpsSetIcntl(F,7,icntl_7));
63     }
64     break;
65 #endif
66   default:
67     if (rank == 0) printf(" PETSC LU:\n");
68     CHKERRQ(MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_LU,&F));
69   }
70 
71   info.fill = 5.0;
72   CHKERRQ(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     CHKERRQ(MatLUFactorNumeric(F,A,&info));
77 
78     /* Test MatSolve() */
79     if (testMatSolve) {
80       CHKERRQ(MatSolve(F,b,x));
81 
82       /* Check the residual */
83       CHKERRQ(MatMult(A,x,u));
84       CHKERRQ(VecAXPY(u,-1.0,b));
85       CHKERRQ(VecNorm(u,NORM_INFINITY,&norm));
86       if (norm > tol) {
87         if (rank == 0) {
88           CHKERRQ(PetscPrintf(PETSC_COMM_SELF,"MatSolve: rel residual %g/%g = %g, LU numfact %d\n",norm,Anorm,norm/Anorm,nfact));
89         }
90       }
91     }
92   }
93 
94   /* Free data structures */
95   CHKERRQ(MatDestroy(&A));
96   CHKERRQ(MatDestroy(&F));
97   CHKERRQ(ISDestroy(&perm));
98   CHKERRQ(ISDestroy(&iperm));
99   CHKERRQ(VecDestroy(&x));
100   CHKERRQ(VecDestroy(&b));
101   CHKERRQ(VecDestroy(&u));
102   ierr = PetscFinalize();
103   return ierr;
104 }
105