[petsc-users] Call KSPConvergedDefault from Fortran

Barry Smith bsmith at mcs.anl.gov
Tue Oct 20 13:25:41 CDT 2015


   It will be very similar to usage from C. You should be able to do something like

   PetscFortranAddress ctx

   call KSPConvergedDefaultCreate(ctx,ierr)

   call KSPMonitorSet(ksp, MyKSPConverged,ctx,kspconvergeddefaultdestroy,ierr)

   then in your MyKSPConverged pass the argument called "dummy" to the default monitor.

  Barry





> On Oct 20, 2015, at 3:36 AM, Natacha BEREUX <natacha.bereux at gmail.com> wrote:
> 
> Dear PETSc users, 
> 
> 
> I am using PETSc  Fortran API. 
> I try to solve a linear system, and would like to define my own convergence test. 
> More precisely, my convergence test calls the default convergence test 
> and then checks the true residual norm : 
>      subroutine MyKSPConverged(ksp,n,rnorm,flag,dummy,ierr)
> 
>       implicit none
> 
> #include <petsc/finclude/petscsys.h>
> #include <petsc/finclude/petscvec.h>
> #include <petsc/finclude/petscksp.h>
> 
>       KSP              ksp
>       PetscErrorCode ierr
>       PetscInt n,dummy
>       KSPConvergedReason flag
>       PetscReal rnorm, true_rnorm
>       Vec true_res
>       !
>       call KSPConvergedDefault(ksp,n,rnorm,flag,PETSC_NULL_OBJECT,ierr)
>       !
>       ! If convergence test succeeds 
>        if ( (flag == KSP_CONVERGED_ATOL) .or.                           &
>      &  (flag == KSP_CONVERGED_RTOL)) then 
>         ! Compute true residual
>         call KSPBuildResidual( ksp , PETSC_NULL_OBJECT,                 &
>      & PETSC_NULL_OBJECT, true_res, ierr )
>         ! Compute true residual norm 
>         call VecNorm(true_res,NORM_2,true_rnorm,ierr)
>         ! And check again convergence with respect to the true residual norm
>         call KSPConvergedDefault(ksp,n,true_rnorm,flag,                 &
>      & PETSC_NULL_OBJECT, ierr)
>         ! Free memory
>         call VecDestroy(true_res, ierr) 
>       endif 
>       !       
>       ierr = 0
> 
>       end
> 
> I get the following error message
> 
> [0]PETSC ERROR: --------------------- Error Message --------------------------------------------------------------
> [0]PETSC ERROR: Null argument, when expecting valid pointer
> [0]PETSC ERROR: Convergence context must have been created with KSPConvergedDefaultCreate()
> [0]PETSC ERROR: See http://www.mcs.anl.gov/petsc/documentation/faq.html for trouble shooting.
> [0]PETSC ERROR: Petsc Release Version 3.6.2, Oct, 02, 2015 
> 
> I understand that I should have  call KSPConvergedDefaultCreate() before calling KSConvergedDefault( ....) 
> 
> Am I right ? 
> 
> If so what is the Fortran calling sequence for  KSPConvergedDefaultCreate ? It is supposed to return a valid pointer and I don't succeed in doing so from Fortran. 
> 
> If you have any idea, I'd greatly appreciate it
> Thank you!
> Natacha 
> 
> 
> PS I have attached my code (it's a slightly modified version of ex2f.F showing  the convergence test code )
> <ex2f.F>



More information about the petsc-users mailing list