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