[petsc-users] Using MUMPS solver with PETSc

Hong hzhang at mcs.anl.gov
Sat Mar 5 11:42:45 CST 2016


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*
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20160305/ffb99319/attachment-0001.html>


More information about the petsc-users mailing list