xref: /petsc/src/mat/tests/ex130.c (revision e00437b97e7e022a766c323f69bd5255763f3ff2)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests external direct solvers. Simplified from ex125.c\n\
3c4762a1bSJed Brown Example: mpiexec -n <np> ./ex130 -f <matrix binary file> -mat_solver_type 1 -mat_superlu_equil \n\n";
4c4762a1bSJed Brown 
5c4762a1bSJed Brown #include <petscmat.h>
6c4762a1bSJed Brown 
7c4762a1bSJed Brown int main(int argc,char **args)
8c4762a1bSJed Brown {
9c4762a1bSJed Brown   Mat            A,F;
10c4762a1bSJed Brown   Vec            u,x,b;
11c4762a1bSJed Brown   PetscMPIInt    rank,size;
12c4762a1bSJed Brown   PetscInt       m,n,nfact,ipack=0;
13c4762a1bSJed Brown   PetscReal      norm,tol=1.e-12,Anorm;
14c4762a1bSJed Brown   IS             perm,iperm;
15c4762a1bSJed Brown   MatFactorInfo  info;
16c4762a1bSJed Brown   PetscBool      flg,testMatSolve=PETSC_TRUE;
17c4762a1bSJed Brown   PetscViewer    fd;              /* viewer */
18c4762a1bSJed Brown   char           file[PETSC_MAX_PATH_LEN]; /* input file name */
19c4762a1bSJed Brown 
209566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
219566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
229566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
23c4762a1bSJed Brown 
24c4762a1bSJed Brown   /* Determine file from which we read the matrix A */
259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&flg));
2628b400f6SJacob Faibussowitsch   PetscCheck(flg,PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f option");
27c4762a1bSJed Brown 
28c4762a1bSJed Brown   /* Load matrix A */
299566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd));
309566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
319566063dSJacob Faibussowitsch   PetscCall(MatLoad(A,fd));
329566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD,&b));
339566063dSJacob Faibussowitsch   PetscCall(VecLoad(b,fd));
349566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&fd));
359566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(A,&m,&n));
36*e00437b9SBarry Smith   PetscCheck(m == n,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ, "This example is not intended for rectangular matrices (%d, %d)", m, n);
379566063dSJacob Faibussowitsch   PetscCall(MatNorm(A,NORM_INFINITY,&Anorm));
38c4762a1bSJed Brown 
39c4762a1bSJed Brown   /* Create vectors */
409566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(b,&x));
419566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&u)); /* save the true solution */
42c4762a1bSJed Brown 
43c4762a1bSJed Brown   /* Test LU Factorization */
449566063dSJacob Faibussowitsch   PetscCall(MatGetOrdering(A,MATORDERINGNATURAL,&perm,&iperm));
45c4762a1bSJed Brown 
469566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-mat_solver_type",&ipack,NULL));
47c4762a1bSJed Brown   switch (ipack) {
48c4762a1bSJed Brown   case 1:
49c4762a1bSJed Brown #if defined(PETSC_HAVE_SUPERLU)
50dd400576SPatrick Sanan     if (rank == 0) printf(" SUPERLU LU:\n");
519566063dSJacob Faibussowitsch     PetscCall(MatGetFactor(A,MATSOLVERSUPERLU,MAT_FACTOR_LU,&F));
52c4762a1bSJed Brown     break;
53c4762a1bSJed Brown #endif
54c4762a1bSJed Brown   case 2:
55c4762a1bSJed Brown #if defined(PETSC_HAVE_MUMPS)
56dd400576SPatrick Sanan     if (rank == 0) printf(" MUMPS LU:\n");
579566063dSJacob Faibussowitsch     PetscCall(MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_LU,&F));
58c4762a1bSJed Brown     {
59c4762a1bSJed Brown       /* test mumps options */
60c4762a1bSJed Brown       PetscInt icntl_7 = 5;
619566063dSJacob Faibussowitsch       PetscCall(MatMumpsSetIcntl(F,7,icntl_7));
62c4762a1bSJed Brown     }
63c4762a1bSJed Brown     break;
64c4762a1bSJed Brown #endif
65c4762a1bSJed Brown   default:
66dd400576SPatrick Sanan     if (rank == 0) printf(" PETSC LU:\n");
679566063dSJacob Faibussowitsch     PetscCall(MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_LU,&F));
68c4762a1bSJed Brown   }
69c4762a1bSJed Brown 
70c4762a1bSJed Brown   info.fill = 5.0;
719566063dSJacob Faibussowitsch   PetscCall(MatLUFactorSymbolic(F,A,perm,iperm,&info));
72c4762a1bSJed Brown 
73c4762a1bSJed Brown   for (nfact = 0; nfact < 1; nfact++) {
74dd400576SPatrick Sanan     if (rank == 0) printf(" %d-the LU numfactorization \n",nfact);
759566063dSJacob Faibussowitsch     PetscCall(MatLUFactorNumeric(F,A,&info));
76c4762a1bSJed Brown 
77c4762a1bSJed Brown     /* Test MatSolve() */
78c4762a1bSJed Brown     if (testMatSolve) {
799566063dSJacob Faibussowitsch       PetscCall(MatSolve(F,b,x));
80c4762a1bSJed Brown 
81c4762a1bSJed Brown       /* Check the residual */
829566063dSJacob Faibussowitsch       PetscCall(MatMult(A,x,u));
839566063dSJacob Faibussowitsch       PetscCall(VecAXPY(u,-1.0,b));
849566063dSJacob Faibussowitsch       PetscCall(VecNorm(u,NORM_INFINITY,&norm));
85c4762a1bSJed Brown       if (norm > tol) {
86dd400576SPatrick Sanan         if (rank == 0) {
879566063dSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_SELF,"MatSolve: rel residual %g/%g = %g, LU numfact %d\n",norm,Anorm,norm/Anorm,nfact));
88c4762a1bSJed Brown         }
89c4762a1bSJed Brown       }
90c4762a1bSJed Brown     }
91c4762a1bSJed Brown   }
92c4762a1bSJed Brown 
93c4762a1bSJed Brown   /* Free data structures */
949566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
959566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&F));
969566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&perm));
979566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&iperm));
989566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
999566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&b));
1009566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&u));
1019566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
102b122ec5aSJacob Faibussowitsch   return 0;
103c4762a1bSJed Brown }
104