[petsc-users] Using MUMPS solver with PETSc

Jose E. Roman jroman at dsic.upv.es
Sat Mar 5 12:01:29 CST 2016


Also, if you call this subroutine several times, make sure to call KSPDestroy() before exit. Alternatively, reuse the KSP object from one call to the next.

Jose


> El 5 mar 2016, a las 18:42, Hong <hzhang at mcs.anl.gov> escribió:
> 
> Ajit :
> I am trying to use PETSc-MUMPS solver to solve linear problem of type "A x = b " 
> I have wrote a subroutine in Fortran 90 for uni-processor. given below.
> I am calling this subroutine inside iterative solver. 
> Code looks fine. 
> 
> This solver is very fast for first few iterations but then becomes slow
> Run your code with runtime option '-log_summary', which displays which routine dominates computation.
>  
> and if matrix is of large-scale, it becomes very-very slow.
> Direct solvers are not scalable in general. Matrix factors might become very dense. You may experiment different matrix orderings. Again, run with '-log_summary'.
> 
> Hong 
> 
> If somebody can help me to understand, if I am doing something wrong here?
> 
> Thanks you. 
> 
> !!!*************************************************
> !!! Subroutine
> SUBROUTINE PETScMUMPS(ksp,pc,A,x,b,ierr)
> 
> #include <petsc/finclude/petscsys.h>
> #include <petsc/finclude/petscvec.h>
> #include <petsc/finclude/petscmat.h>
> #include <petsc/finclude/petscksp.h>
> #include <petsc/finclude/petscpc.h>
> 
> Vec              x,b
> Mat              A
> KSP              ksp
> PC               pc
> 
> Mat              F
> PetscInt         ival,icntl
> 
> PetscErrorCode   ierr
> 
>     call KSPCreate(PETSC_COMM_SELF,ksp,ierr)
>     call KSPSetOperators(ksp,A,A,ierr)
>     call KSPSetType(ksp,KSPPREONLY,ierr)
>     call KSPGetPC(ksp,pc,ierr)
> 
>     !call PCSetType(pc,PCLU,ierr)           !! LU Factorization
>     call PCSetType(pc,PCCHOLESKY,ierr)     !! Cholesky Factorization
> 
>     call PCFactorSetMatSolverPackage(pc,MATSOLVERMUMPS,ierr)
>     call PCFactorSetUpMatSolverPackage(pc,ierr)
>     call PCFactorGetMatrix(pc,F,ierr)
> 
>     !! sequential ordering
>     icntl = 7
>     ival  = 2
>     call MatMumpsSetIcntl(F,icntl,ival,ierr)
> 
>     call KSPSetFromOptions(ksp,ierr)
>     call KSPGetPC(ksp,pc,ierr)
> 
>     call KSPSolve(ksp,x,b,ierr)
> 
> END SUBROUTINE PETScMUMPS
> !!!*************************************************
> 
> 
> Ajit Desai
> --
>   PhD Scholar, Carleton University 
>   Ottawa, Canada
> 



More information about the petsc-users mailing list