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