call PetscOptionsSetValue in fortran

Barry Smith bsmith at mcs.anl.gov
Fri Jul 17 10:51:53 CDT 2009


    You are calling KSPSetFromOptions() BEFORE setting the PC type to  
hypre and setting boomeramg and before you call PetscOptionsSetValue().
You should call PetscOptionsSetValue() then PCSetType() then  
PCHypreSetType() then KSPSetFromOptions().

    Barry

On Jul 17, 2009, at 4:52 AM, Klaij, Christiaan wrote:

>
> I'm trying to change the options of the Hypre preconditioner using  
> PetscOptionsSetValue in a fortran program, but I must be doing  
> something wrong, see the session below. It works fine from the  
> command line, though. As an example, I took ex12f from src/ksp/ksp/ 
> examples/tests (petsc-2.3.3-p13) and modified it a little.
>
>
> $ cat ex12f.F
> !
>       program main
>        implicit none
>
> #include "include/finclude/petsc.h"
> #include "include/finclude/petscvec.h"
> #include "include/finclude/petscmat.h"
> #include "include/finclude/petscpc.h"
> #include "include/finclude/petscksp.h"
> #include "include/finclude/petscviewer.h"
> !
> !  This example is the Fortran version of ex6.c.  The program reads  
> a PETSc matrix
> !  and vector from a file and solves a linear system.  Input  
> arguments are:
> !        -f <input_file> : file to load.  For a 5X5 example of the 5- 
> pt. stencil
> !                          use the file petsc/src/mat/examples/ 
> matbinary.ex
> !
>
>       PetscErrorCode  ierr
>       PetscInt its
>       PetscTruth flg
>       PetscScalar      norm,none
>       Vec              x,b,u
>       Mat              A
>       character*(128)  f
>       PetscViewer      fd
>       MatInfo          info(MAT_INFO_SIZE)
>       KSP              ksp
> ! cklaij: adding pc
>       PC               pc
> ! cklaij: adding pc end
>
>       none = -1.0
>       call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
>
> ! Read in matrix and RHS
>       call PetscOptionsGetString(PETSC_NULL_CHARACTER,'-f',f,flg,ierr)
>       print *,f
>       call  
> PetscViewerBinaryOpen(PETSC_COMM_WORLD,f,FILE_MODE_READ,     &
>      &     fd,ierr)
>
>       call MatLoad(fd,MATSEQAIJ,A,ierr)
>
> ! Get information about matrix
>       call MatGetInfo(A,MAT_GLOBAL_SUM,info,ierr)
>       write(*,100)  
> info(MAT_INFO_ROWS_GLOBAL),                          &
>      &   
> info(MAT_INFO_COLUMNS_GLOBAL),                                  &
>      &   
> info(MAT_INFO_ROWS_LOCAL),info(MAT_INFO_COLUMNS_LOCAL),         &
>      &   
> info(MAT_INFO_BLOCK_SIZE),info(MAT_INFO_NZ_ALLOCATED),          &
>      &   
> info(MAT_INFO_NZ_USED),info(MAT_INFO_NZ_UNNEEDED),              &
>      &   
> info(MAT_INFO_MEMORY),info(MAT_INFO_ASSEMBLIES),                &
>      &  info(MAT_INFO_MALLOCS)
>
>  100  format(11(g7.1,1x))
>       call VecLoad(fd,PETSC_NULL_CHARACTER,b,ierr)
>       call PetscViewerDestroy(fd,ierr)
>
> ! Set up solution
>       call VecDuplicate(b,x,ierr)
>       call VecDuplicate(b,u,ierr)
>
> ! Solve system
>       call KSPCreate(PETSC_COMM_WORLD,ksp,ierr)
>       call KSPSetOperators(ksp,A,A,DIFFERENT_NONZERO_PATTERN,ierr)
>       call KSPSetFromOptions(ksp,ierr)
> ! cklaij: try boomeramg
>       call KSPGetPC(ksp,pc,ierr)
>       call PCSetType(pc,PCHYPRE,ierr)
>       call PCHYPRESetType(pc,"boomeramg",ierr)
>       call PetscOptionsSetValue
>      &     ("-pc_hypre_boomeramg_strong_threshold","0.5",ierr)
> ! cklaij: try boomeramg end
>       call KSPSolve(ksp,b,x,ierr)
>
> ! Show result
>       call MatMult(A,x,u,ierr)
>       call VecAXPY(u,none,b,ierr)
>       call VecNorm(u,NORM_2,norm,ierr)
>       call KSPGetIterationNumber(ksp,its,ierr)
>       print*, 'Number of iterations = ',its
>       print*, 'Residual norm = ',norm
>
> ! Cleanup
>       call KSPDestroy(ksp,ierr)
>       call VecDestroy(b,ierr)
>       call VecDestroy(x,ierr)
>       call VecDestroy(u,ierr)
>       call MatDestroy(A,ierr)
>
>       call PetscFinalize(ierr)
>       end
>
> $ ./ex12f  -f ../../../../mat/examples/matbinary.ex -ksp_view | grep  
> Threshold
>     HYPRE BoomerAMG: Threshold for strong coupling 0.25
> $ ./ex12f  -f ../../../../mat/examples/matbinary.ex -pc_type hypre - 
> pc_hypre_type boomeramg -pc_hypre_boomeramg_strong_threshold 0.5 - 
> ksp_view | grep Threshold
>     HYPRE BoomerAMG: Threshold for strong coupling 0.5
> $
>
>
> <mime-attachment.jpeg><mime-attachment.jpeg>
> dr. ir. Christiaan Klaij
> CFD Researcher
> Research & Development
> MARIN
> 2, Haagsteeg
> c.klaij at marin.nl
> P.O. Box 28
> T +31 317 49 39 11
> 6700 AA  Wageningen
> F +31 317 49 32 45
> T  +31 317 49 33 44
> The Netherlands
> I  www.marin.nl
>
>
> MARIN webnews: First AMT'09 conference, Nantes, France, September 1-2
>
>
> This e-mail may be confidential, privileged and/or protected by  
> copyright. If you are not the intended recipient, you should return  
> it to the sender immediately and delete your copy from your system.
>



More information about the petsc-users mailing list