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 ierr = MPI_Comm_rank(PETSC_COMM_WORLD, &rank);CHKERRMPI(ierr); 23 ierr = MPI_Comm_size(PETSC_COMM_WORLD, &size);CHKERRMPI(ierr); 24 25 /* Determine file from which we read the matrix A */ 26 ierr = PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&flg);CHKERRQ(ierr); 27 if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER,"Must indicate binary file with the -f option"); 28 29 /* Load matrix A */ 30 ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd);CHKERRQ(ierr); 31 ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 32 ierr = MatLoad(A,fd);CHKERRQ(ierr); 33 ierr = VecCreate(PETSC_COMM_WORLD,&b);CHKERRQ(ierr); 34 ierr = VecLoad(b,fd);CHKERRQ(ierr); 35 ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr); 36 ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr); 37 if (m != n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ, "This example is not intended for rectangular matrices (%d, %d)", m, n); 38 ierr = MatNorm(A,NORM_INFINITY,&Anorm);CHKERRQ(ierr); 39 40 /* Create vectors */ 41 ierr = VecDuplicate(b,&x);CHKERRQ(ierr); 42 ierr = VecDuplicate(x,&u);CHKERRQ(ierr); /* save the true solution */ 43 44 /* Test LU Factorization */ 45 ierr = MatGetOrdering(A,MATORDERINGNATURAL,&perm,&iperm);CHKERRQ(ierr); 46 47 ierr = PetscOptionsGetInt(NULL,NULL,"-mat_solver_type",&ipack,NULL);CHKERRQ(ierr); 48 switch (ipack) { 49 case 1: 50 #if defined(PETSC_HAVE_SUPERLU) 51 if (!rank) printf(" SUPERLU LU:\n"); 52 ierr = MatGetFactor(A,MATSOLVERSUPERLU,MAT_FACTOR_LU,&F);CHKERRQ(ierr); 53 break; 54 #endif 55 case 2: 56 #if defined(PETSC_HAVE_MUMPS) 57 if (!rank) printf(" MUMPS LU:\n"); 58 ierr = MatGetFactor(A,MATSOLVERMUMPS,MAT_FACTOR_LU,&F);CHKERRQ(ierr); 59 { 60 /* test mumps options */ 61 PetscInt icntl_7 = 5; 62 ierr = MatMumpsSetIcntl(F,7,icntl_7);CHKERRQ(ierr); 63 } 64 break; 65 #endif 66 default: 67 if (!rank) printf(" PETSC LU:\n"); 68 ierr = MatGetFactor(A,MATSOLVERPETSC,MAT_FACTOR_LU,&F);CHKERRQ(ierr); 69 } 70 71 info.fill = 5.0; 72 ierr = MatLUFactorSymbolic(F,A,perm,iperm,&info);CHKERRQ(ierr); 73 74 for (nfact = 0; nfact < 1; nfact++) { 75 if (!rank) printf(" %d-the LU numfactorization \n",nfact); 76 ierr = MatLUFactorNumeric(F,A,&info);CHKERRQ(ierr); 77 78 /* Test MatSolve() */ 79 if (testMatSolve) { 80 ierr = MatSolve(F,b,x);CHKERRQ(ierr); 81 82 /* Check the residual */ 83 ierr = MatMult(A,x,u);CHKERRQ(ierr); 84 ierr = VecAXPY(u,-1.0,b);CHKERRQ(ierr); 85 ierr = VecNorm(u,NORM_INFINITY,&norm);CHKERRQ(ierr); 86 if (norm > tol) { 87 if (!rank) { 88 ierr = PetscPrintf(PETSC_COMM_SELF,"MatSolve: rel residual %g/%g = %g, LU numfact %d\n",norm,Anorm,norm/Anorm,nfact);CHKERRQ(ierr); 89 } 90 } 91 } 92 } 93 94 /* Free data structures */ 95 ierr = MatDestroy(&A);CHKERRQ(ierr); 96 ierr = MatDestroy(&F);CHKERRQ(ierr); 97 ierr = ISDestroy(&perm);CHKERRQ(ierr); 98 ierr = ISDestroy(&iperm);CHKERRQ(ierr); 99 ierr = VecDestroy(&x);CHKERRQ(ierr); 100 ierr = VecDestroy(&b);CHKERRQ(ierr); 101 ierr = VecDestroy(&u);CHKERRQ(ierr); 102 ierr = PetscFinalize(); 103 return ierr; 104 } 105