[petsc-users] Issues with subsequent GMRES Solves
C.A.(Sandy) Mader
mader at utias.utoronto.ca
Thu May 13 09:20:49 CDT 2010
I was able to get a version of the code working last night by destroying
all of the PETSc objects between subsequent adjoint setup calls. Based
on some of my previous results I suspect that I am having an issue with
resetting the matrices used in the linear solve. I am running some more
tests now to confirm this...
Sandy
On Wed, 2010-05-12 at 14:02 -0500, Barry Smith wrote:
> On May 12, 2010, at 1:55 PM, C.A.(Sandy) Mader wrote:
>
> > See responses in text below:
> > On Wed, 2010-05-12 at 13:05 -0500, Barry Smith wrote:
> >> Are you calling KSPSetOperators(ksp,A,B,DIFFERENT_NONZERO_PATTERN); before the second adjoint solve? Or are you trying to reuse a preconditioner from the first set of solves?
> > I am calling:
> > call KSPSetOperators(ksp,dRdWT,dRdWPreT, &
> > DIFFERENT_NONZERO_PATTERN,PETScIerr)
>
> When you use this flag it totally rebuilds the preconditioner, there is nothing else you need to do, so there is no "reason" that it shouldn't behave the same whether it is the first or second time you do a solve. Note: you can always destroy the KSP and make a new one to see if that makes a difference (it should not).
>
> Have you run the code with valgrind http://www.mcs.anl.gov/petsc/petsc-as/documentation/faq.html#valgrind to make sure there is no subtle memory corruption causing problems.
>
> Is the first residual norm of the "good" and bad second adjoint iteration the same? What about the second?
>
> You can use MatView() with a binary viewer to save all the matrices for the "good" and bad second adjoint, load them into Matlab and see if they are the same. Same with the right hand side.
>
>
> Barry
>
> >>
> >> When the "second adjoint" converges and when it doesn't converge are you using the exact same matrix A, matrix B and right hand side?
> > As far as I know yes. For the small case (the one that works), I printed
> > them out and compared them and they were identical. I haven't been able
> > to run the same comparison for the larger case...
> >>
> >> Barry
> >>
> >>
> >>
> >> On May 12, 2010, at 1:00 PM, C.A.(Sandy) Mader wrote:
> >>
> >>> See responses in text below:
> >>> On Wed, 2010-05-12 at 12:48 -0500, Barry Smith wrote:
> >>>> Is it the same matrix or a different matrix that is solved with the "second adjoint system"? Or is it just a different right hand side?
> >>> It is a different matrix. In this case I solve the CFD at a given flow
> >>> condition, then solve the adjoint twice (two different RHS) both
> >>> converge. Then I solve the CFD case for a second flow condition,
> >>> followed by two more adjoint solutions. Neither of which solve. However,
> >>> if I solver directly for the second flow condition, the adojint system
> >>> will solve...
> >>>>
> >>>> If it is a different system how do you know it will even converge for the new matrix? Maybe the new matrix is much more difficult?
> >>> see above answer, if the same system will solve under one set of
> >>> conditions but not under another...
> >>>>
> >>>> Is it possible the "second adjoint system" is singular?
> >>> The second system will solve if I solve it on its own.
> >>>
> >>>> Given that this is a very strong preconditioner the convergence for the first adjoint system is pretty poor.
> >>> We use an approximate preconditioner to reduce the memory requirements
> >>> of the solver. I do not have this problem if I use the exact matrix as
> >>> the preconditioner but it limits my problem size....
> >>>>
> >>>> Does this happen also for much smaller problems? Recommend running much smaller problem to see what happens convergence wise, then run that much smaller problem with -pc_type lu to see if the direct solver handles it happily.
> >>> I have tried it on a smaller case and this problem does not occur...
> >>>>
> >>>>
> >>>>
> >>>> On May 12, 2010, at 12:39 PM, C.A.(Sandy) Mader wrote:
> >>>>
> >>>>> I am using PETSc to solve an adjoint system for a CFD code.
> >>>>> The issue I am having is that if I try to run successive cfd/adjoint
> >>>>> solutions, the adjoint fails to converge after the first cfd/adjoint
> >>>>> solution cycle. i.e.
> >>>>> Solve CFD (ok)
> >>>>> Solve adjoint system with multiple RHS's (ok)
> >>>>> Solve CFD at new point (ok)
> >>>>> Solve adjoint system at new point with multiple RHS's (GMRES iterations
> >>>>> stall)
> >>>>>
> >>>>> PETSc is not used to solve the CFD system, only the adjoint system.
> >>>>>
> >>>>> I am currently using the following solver options:
> >>>>> 'Adjoint solver type': 'GMRES',\
> >>>>> 'adjoint relative tolerance':1e-6,\
> >>>>> 'adjoint absolute tolerance':1e-16,\
> >>>>> 'adjoint max iterations': 1500,\
> >>>>> 'adjoint restart iteration' : 80,\
> >>>>> 'adjoint monitor step': 10,\
> >>>>> 'Preconditioner Side': 'right',\
> >>>>> 'Matrix Ordering': 'NestedDissection',\
> >>>>> 'Global Preconditioner Type': 'Additive Schwartz',\
> >>>>> 'Local Preconditioner Type' : 'ILU',\
> >>>>> 'ILU Fill Levels': 3,\
> >>>>> 'ASM Overlap' : 5,\
> >>>>>
> >>>>> Also, I am using the fortran interface to PETSc.
> >>>>>
> >>>>> As a sample, I have included two convergence histories. Both are at the
> >>>>> same converged CFD point(this case is ~1.2million cell so the PETSc
> >>>>> system is approx 6million degrees of freedom). The first is an example
> >>>>> where the adjoint system is solved the first time through the cycle. The
> >>>>> second is an example of where the adjoint is solved the second time
> >>>>> through the cycle:
> >>>>> 1) as first point in cycle
> >>>>> # ... KSP properties:
> >>>>> # type : gmres
> >>>>> # tolerances : rel = 1.0E-10
> >>>>> # abs = 1.0E-16
> >>>>> # div = 1.0E+05
> >>>>> # max.iter. : 1500
> >>>>> # precond.type: asm
> >>>>> Solving ADjoint Transpose with PETSc...
> >>>>> 0 KSP Residual norm 0.1392E+00
> >>>>> 10 KSP Residual norm 0.6394E-01
> >>>>> 20 KSP Residual norm 0.6106E-01
> >>>>> 30 KSP Residual norm 0.6019E-01
> >>>>> 40 KSP Residual norm 0.5941E-01
> >>>>> 50 KSP Residual norm 0.5876E-01
> >>>>> 60 KSP Residual norm 0.5602E-01
> >>>>> 70 KSP Residual norm 0.4915E-01
> >>>>> 80 KSP Residual norm 0.3994E-01
> >>>>> 90 KSP Residual norm 0.3892E-01
> >>>>> 100 KSP Residual norm 0.3854E-01
> >>>>> 110 KSP Residual norm 0.3794E-01
> >>>>> 120 KSP Residual norm 0.3717E-01
> >>>>> 130 KSP Residual norm 0.3630E-01
> >>>>> 140 KSP Residual norm 0.3415E-01
> >>>>> ...
> >>>>> 900 KSP Residual norm 0.2437E-09
> >>>>> 910 KSP Residual norm 0.1452E-09
> >>>>> 920 KSP Residual norm 0.1025E-09
> >>>>> 930 KSP Residual norm 0.6875E-10
> >>>>> 940 KSP Residual norm 0.4141E-10
> >>>>> 950 KSP Residual norm 0.2317E-10
> >>>>> 960 KSP Residual norm 0.1559E-10
> >>>>> Solving ADjoint Transpose with PETSc time (s) = 391.43
> >>>>> Norm of error = 0.1366E-10 Iterations = 965
> >>>>> ------------------------------------------------
> >>>>> PETSc solver converged after 965 iterations.
> >>>>> ------------------------------------------------
> >>>>>
> >>>>> 2) As second point in cycle
> >>>>> # ... KSP properties:
> >>>>> # type : gmres
> >>>>> # tolerances : rel = 1.0E-10
> >>>>> # abs = 1.0E-16
> >>>>> # div = 1.0E+05
> >>>>> # max.iter. : 1500
> >>>>> # precond.type: asm
> >>>>> Solving ADjoint Transpose with PETSc...
> >>>>> 0 KSP Residual norm 0.1392E+00
> >>>>> 10 KSP Residual norm 0.6400E-01
> >>>>> 20 KSP Residual norm 0.6140E-01
> >>>>> 30 KSP Residual norm 0.6060E-01
> >>>>> 40 KSP Residual norm 0.5995E-01
> >>>>> 50 KSP Residual norm 0.5974E-01
> >>>>> 60 KSP Residual norm 0.5971E-01
> >>>>> 70 KSP Residual norm 0.5957E-01
> >>>>> 80 KSP Residual norm 0.5906E-01
> >>>>> 90 KSP Residual norm 0.5906E-01
> >>>>> 100 KSP Residual norm 0.5906E-01
> >>>>> 110 KSP Residual norm 0.5903E-01
> >>>>> 120 KSP Residual norm 0.5901E-01
> >>>>> 130 KSP Residual norm 0.5901E-01
> >>>>> 140 KSP Residual norm 0.5901E-01
> >>>>> ...
> >>>>> 1400 KSP Residual norm 0.5895E-01
> >>>>> 1410 KSP Residual norm 0.5895E-01
> >>>>> 1420 KSP Residual norm 0.5895E-01
> >>>>> 1430 KSP Residual norm 0.5895E-01
> >>>>> 1440 KSP Residual norm 0.5895E-01
> >>>>> 1450 KSP Residual norm 0.5895E-01
> >>>>> 1460 KSP Residual norm 0.5895E-01
> >>>>> 1470 KSP Residual norm 0.5895E-01
> >>>>> 1480 KSP Residual norm 0.5895E-01
> >>>>> 1490 KSP Residual norm 0.5895E-01
> >>>>> 1500 KSP Residual norm 0.5895E-01
> >>>>> Solving ADjoint Transpose with PETSc time (s) = 516.59
> >>>>> Norm of error = 0.5895E-01 Iterations = 1500
> >>>>> ------------------------------------------------
> >>>>> PETSc solver converged after 1500 iterations.
> >>>>> ------------------------------------------------
> >>>>> I have tried both resetting the operators between cycles and completely
> >>>>> destroying the KSP object between cycles. Both give the same result. Is
> >>>>> there a step I am missing to properly reset the system for the
> >>>>> subsequent solves?
> >>>>>
> >>>>> Thanks in advance for your input...
> >>>>>
> >>>>> Sandy
> >>>>> --
> >>>>> C.A.(Sandy) Mader
> >>>>> Ph.D Candidate
> >>>>> Multidisciplinary Design Optimization Laboratory
> >>>>> University of Toronto Institute for Aerospace Studies
> >>>>> mader at utias.utoronto.ca
> >>>>>
> >>>>
> >>>>
> >>>
> >>
> >>
> >
>
>
More information about the petsc-users
mailing list