[petsc-users] monitoring the convergence of fieldsplit 0 and 1

Klaij, Christiaan C.Klaij at marin.nl
Fri Jan 13 03:46:35 CST 2017


Barry,

It's been a while but I'm finally using this function in
3.7.4. Is it supposed to work with fieldsplit? Here's why.

I'm solving a Navier-Stokes system with fieldsplit (pc has one
velocity solve and one pressure solve) and trying to retrieve the
totals like this:

  CALL KSPSolve(ksp_system,rr_system,xx_system,ierr); CHKERRQ(ierr)
  CALL PCFieldSplitGetSubKSP(pc_system,numsplit,subksp,ierr); CHKERRQ(ierr)
  CALL KSPGetTotalIterations(subksp(1),nusediter_vv,ierr); CHKERRQ(ierr)
  CALL KSPGetTotalIterations(subksp(2),nusediter_pp,ierr); CHKERRQ(ierr)
print *, 'nusediter_vv', nusediter_vv
print *, 'nusediter_pp', nusediter_pp

Running the code shows this surprise:

  Linear sys_fieldsplit_0_ solve converged due to CONVERGED_RTOL iterations 1
  Linear sys_fieldsplit_0_ solve converged due to CONVERGED_RTOL iterations 1
  Linear sys_fieldsplit_0_ solve converged due to CONVERGED_RTOL iterations 2
  Linear sys_fieldsplit_0_ solve converged due to CONVERGED_RTOL iterations 2
  Linear sys_fieldsplit_0_ solve converged due to CONVERGED_RTOL iterations 7
  Linear sys_fieldsplit_0_ solve converged due to CONVERGED_RTOL iterations 7
  Linear sys_fieldsplit_0_ solve converged due to CONVERGED_RTOL iterations 7
  Linear sys_fieldsplit_0_ solve converged due to CONVERGED_RTOL iterations 7
  Linear sys_fieldsplit_0_ solve converged due to CONVERGED_RTOL iterations 8
  Linear sys_fieldsplit_0_ solve converged due to CONVERGED_RTOL iterations 8
  Linear sys_fieldsplit_0_ solve converged due to CONVERGED_RTOL iterations 8
  Linear sys_fieldsplit_0_ solve converged due to CONVERGED_RTOL iterations 8

  Linear sys_fieldsplit_1_ solve converged due to CONVERGED_RTOL iterations 22
  Linear sys_fieldsplit_1_ solve converged due to CONVERGED_RTOL iterations 6
  Linear sys_fieldsplit_1_ solve converged due to CONVERGED_RTOL iterations 3
  Linear sys_fieldsplit_1_ solve converged due to CONVERGED_RTOL iterations 2
  Linear sys_fieldsplit_1_ solve converged due to CONVERGED_RTOL iterations 2
  Linear sys_fieldsplit_1_ solve converged due to CONVERGED_RTOL iterations 2

 nusediter_vv          37
 nusediter_pp          37

So the value of nusediter_pp is indeed 37, but for nusediter_vv
it should be 66. Any idea what went wrong?

Chris



dr. ir. Christiaan Klaij  | CFD Researcher | Research & Development
MARIN | T +31 317 49 33 44 | mailto:C.Klaij at marin.nl | http://www.marin.nl

MARIN news: http://www.marin.nl/web/News/News-items/MARIN-wishes-you-a-challenging-inspiring-2017.htm

________________________________________
From: Barry Smith <bsmith at mcs.anl.gov>
Sent: Saturday, April 11, 2015 12:27 AM
To: Klaij, Christiaan
Cc: petsc-users at mcs.anl.gov
Subject: Re: [petsc-users] monitoring the convergence of fieldsplit 0 and 1

  Chris,

  I have added KSPGetTotalIterations() to the branch barry/add-ksp-total-iterations/master and next. After tests it will go into master

  Barry

> On Apr 10, 2015, at 8:07 AM, Klaij, Christiaan <C.Klaij at marin.nl> wrote:
>
> Barry,
>
> Sure, I can call PCFieldSplitGetSubKSP() to get the fieldsplit_0
> ksp and then KSPGetIterationNumber, but what does this number
> mean?
>
> It appears to be the number of iterations of the last time that
> the subsystem was solved, right? If so, this corresponds to the
> last iteration of the coupled system, how about all the previous
> iterations?
>
> Chris
> ________________________________________
> From: Barry Smith <bsmith at mcs.anl.gov>
> Sent: Friday, April 10, 2015 2:48 PM
> To: Klaij, Christiaan
> Cc: petsc-users at mcs.anl.gov
> Subject: Re: [petsc-users] monitoring the convergence of fieldsplit 0 and 1
>
>   Chris,
>
>     It appears you should call PCFieldSplitGetSubKSP() and then get the information you want out of the individual KSPs. If this doesn't work please let us know.
>
>   Barry
>
>> On Apr 10, 2015, at 6:48 AM, Klaij, Christiaan <C.Klaij at marin.nl> wrote:
>>
>> A question when using PCFieldSplit: for each linear iteration of
>> the system, how many iterations for fielsplit 0 and 1?
>>
>> One way to find out is to run with -ksp_monitor,
>> -fieldsplit_0_ksp_monitor and -fieldsplit_0_ksp_monitor. This
>> gives the complete convergence history.
>>
>> Another way, suggested by Matt, is to use -ksp_monitor,
>> -fieldsplit_0_ksp_converged_reason and
>> -fieldsplit_1_ksp_converged_reason. This gives only the totals
>> for fieldsplit 0 and 1 (but without saying for which one).
>>
>> Both ways require to somehow process the output, which is a bit
>> inconvenient. Could KSPGetResidualHistory perhaps return (some)
>> information on the subsystems' convergence for processing inside
>> the code?
>>
>> Chris
>>
>>
>> dr. ir. Christiaan Klaij
>> CFD Researcher
>> Research & Development
>> E mailto:C.Klaij at marin.nl
>> T +31 317 49 33 44
>>
>>
>> MARIN
>> 2, Haagsteeg, P.O. Box 28, 6700 AA Wageningen, The Netherlands
>> T +31 317 49 39 11, F +31 317 49 32 45, I www.marin.nl
>>
>



More information about the petsc-users mailing list