[petsc-users] Help with KSPSetConvergenceTest

Barry Smith bsmith at petsc.dev
Sun May 7 10:47:07 CDT 2023


  Note, you must call the default test on iteration 0. This is how it determines the initial residual for relative residual tests etc. 

  I recommend not having the if () test at all. Instead, always call the default convergence test first and then change the flag value it provides if needed before doing any custom checking.

> On May 7, 2023, at 11:38 AM, Barry Smith <bsmith at petsc.dev> wrote:
> 
> 
>    The code will not work as written because KSPConvergedDefault() requires a context created with KSPConvergedDefaultCreate(). 
> 
>    Here is a starting point for what you need, in main
> 
>    integer*8 defaultctx
>    extern MyKSPConverged, KSPConvergedDefaultDestroy
> 
> KSPDefaultConvergedCreate(defaultctx,ierr)
> KSPSetConvergenceTest <https://petsc.org/main/manualpages/KSP/KSPSetConvergenceTest/>(ksp,MyKSPConverged,defaultctx,KSPConvergedDefaultDestroy, ierr)
> 
> 
> subroutine MyKSPConverged(ksp,n,rnorm,flag,defaultctx,ierr)
> 
>        KSP              ksp
>        PetscErrorCode ierr
>        PetscInt n
>        integer*8 defaultctx
>        KSPConvergedReason flag
>        PetscReal rnorm
> 
>        if (n>1) then
>          call KSPConvergedDefault(ksp, n, rnorm, flag,defaultctx, ierr)
>        else
>          flag = 0
>        endif
>        ierr = 0
> 
>     end subroutine MyKSPConverged
> 
> 
> 
>> On May 7, 2023, at 10:09 AM, Matthew Knepley <knepley at gmail.com> wrote:
>> 
>> On Sun, May 7, 2023 at 10:02 AM Edoardo alinovi <edoardo.alinovi at gmail.com <mailto:edoardo.alinovi at gmail.com>> wrote:
>>> Thanks, 
>>> 
>>> Is this a reasonable thing to do if I want to replicate what KSP is doing by default?
>> 
>> Yes. The other option is to pass along 'dummy'
>> 
>>   Thanks,
>> 
>>     Matt
>> 
>> -- 
>> What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.
>> -- Norbert Wiener
>> 
>> https://www.cse.buffalo.edu/~knepley/ <http://www.cse.buffalo.edu/~knepley/>
> 

-------------- next part --------------
An HTML attachment was scrubbed...
URL: <http://lists.mcs.anl.gov/pipermail/petsc-users/attachments/20230507/4a885e47/attachment-0003.html>


More information about the petsc-users mailing list