[petsc-users] Reuse symbolic factorization with petsc - mumps

Thibault Faney tibo at berkeley.edu
Thu Mar 14 14:52:16 CDT 2013


Hong:

thank you for your answer.
I can replace DIFFERENT_NONZERO_PATTERN by SAME_NONZERO_PATTERN, however
this will give me the wrong answer as the zero pattern changes from a time
step to another (e.g. the matrix A at time step 1 has more zeros than at
time step 2) . What I can precompute is the maximal non-zero pattern. Do
you mean I should compute the facorization for the maximal zero pattern
before time iterations and then use SAME_NONZERO_PATTERN during the
iterations ?


On Thu, Mar 14, 2013 at 3:29 PM, Hong Zhang <hzhang at mcs.anl.gov> wrote:

> Tibo:
>
> Replace
> call KSPSetOperators(KSP_A, A, A, DIFFERENT_NONZERO_PATTERN, ierr)
> with
> call KSPSetOperators(KSP_A, A, A, SAME_NONZERO_PATTERN, ierr)
>
> Run your code with option '-log_summary' to see the number of calling
> MatLUSymbolic() and MatLUNumeric().
>
> Hong
>
> >
> > I am using the cray-petsc implementation on the NERSC system.
> Specifically,
> > I have implemented a runge kutta solver, and I need to solve a linear
> system
> > Ax = b several times per time step. I am using PETSC to interface with
> the
> > MUMPS solver for direct factorization.
> >
> > The matrix A is constant over a whole time step, but changes at each time
> > step. What I do now is that at the beginning of each time step, I compute
> > the matrix A, use the MUMPS solver to compute a LU factorization of A,
> and
> > then use KSPSolve to solve each linear system during the same time step.
> > Here is an example of the succession of calls (inspired from example
> > ex52.c):
> >
> > Once per time step:
> > call KSPSetOperators(KSP_A, A, A, DIFFERENT_NONZERO_PATTERN, ierr)
> > call KSPSetType(KSP_A, KSPPREONLY, ierr);
> > call KSPGetPC(KSP_A, PC_A, ierr)
> > call KSPSetTolerances(KSP_A, 1.d-20, PETSC_DEFAULT_DOUBLE_PRECISION,
> > PETSC_DEFAULT_DOUBLE_PRECISION, PETSC_DEFAULT_INTEGER, ierr)
> > call PCSetType(PC_A, PCLU, ierr)
> > call PCFactorSetMatSolverPackage(PC_A, MATSOLVERMUMPS, ierr)
> > call PCFactorSetUpMatSolverPackage(PC_A, ierr)
> > call PCFactorGetMatrix(PC_A, PC_A, ierr)
> > call MatMumpsSetIcntl(PC_A, 7, 5, ierr)
> > call PCSetup(PC_A, ierr)
> > call KSPSetUp(KSP_A, ierr)
> >
> > Then several times per time step:
> > call KSPSolve( KSP_A, b, x, ierr)
> >
> > This works fine.
> > However, I am able to precompute the maximal sparsity pattern of A before
> > the time iterations. Therefore, I would like to use this information by
> > precomputing the symbolic factorization before the time iterations start,
> > and only compute the numerical factorization at each time step. Is there
> a
> > way to use MUMPS in PETSC to do this ? I would appreciate any
> help/reference
> > on this topic.
> >
> > Thank you,
> >
> > Tibo
>
-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20130314/fabc8fa3/attachment.html>


More information about the petsc-users mailing list