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