<div dir="ltr"><div class="gmail_quote"><div dir="ltr">On Tue, Sep 4, 2018 at 10:12 AM Kenneth Hall, SC.D. <<a href="mailto:kenneth.c.hall@duke.edu">kenneth.c.hall@duke.edu</a>> wrote:<br></div><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex">Hi All,<br>
<br>
I am solving a nonlinear problem using a matrix free approach with Petsc. My code is written in Fortran90.<br>
I have read the documentation, and have had some success, but I am having problems with two features.<br>
<br>
1. The code runs fine with no preconditioning.  But when I using the -user_precond option, my preconditioner never gets called by Petsc, even though I make the appropriate calls to PCSetType and PCShellSetApply.<br></blockquote><div><br></div><div>Lets do this first. I think your PC choice is getting overwritten by options. Can you run with -snes_view and</div><div>send the output when trying to use your MF PC?</div><div><br></div><div>  Thanks,</div><div><br></div><div>     Matt</div><div> </div><blockquote class="gmail_quote" style="margin:0 0 0 .8ex;border-left:1px #ccc solid;padding-left:1ex">
2. I am using the GMRES solver, which defaults to the L2-norm minimization line search (SNESLINESEARCHL2). I would like to use backtracking (SNESLINESEARCHBT).  I have specified backtracking using SNESLineSearchSetType, but I then get the following error:<br>
<br>
[0]PETSC ERROR: Object is in wrong state<br>
[0]PETSC ERROR: MatMFFDSetBase() has not been called, this is often caused by forgetting to call <br>
                MatAssemblyBegin/End on the first Mat in the SNES compute function<br>
[0]PETSC ERROR: See <a href="http://www.mcs.anl.gov/petsc/documentation/faq.html" rel="noreferrer" target="_blank">http://www.mcs.anl.gov/petsc/documentation/faq.html</a> for trouble shooting.<br>
[0]PETSC ERROR: Petsc Release Version 3.9.3, Jul, 02, 2018 <br>
[0]PETSC ERROR: mustang on a arch-darwin-c-debug-double-real named Macintosh-11.local by hall Tue Sep  4 07:29:25 2018<br>
[0]PETSC ERROR: Configure options --with-scalar-type=real<br>
[0]PETSC ERROR: #1 MatMult_MFFD() line 314 in /Users/hall/Documents/Fortran_Codes/Packages/petsc-3.9.3/src/mat/impls/mffd/mffd.c<br>
[0]PETSC ERROR: #2 MatMult() line 2311 in /Users/hall/Documents/Fortran_Codes/Packages/petsc-3.9.3/src/mat/interface/matrix.c<br>
[0]PETSC ERROR: #3 SNESLineSearchApply_BT() line 130 in /Users/hall/Documents/Fortran_Codes/Packages/petsc-3.9.3/src/snes/linesearch/impls/bt/linesearchbt.c<br>
[0]PETSC ERROR: #4 SNESLineSearchApply() line 645 in /Users/hall/Documents/Fortran_Codes/Packages/petsc-3.9.3/src/snes/linesearch/interface/linesearch.c<br>
[0]PETSC ERROR: #5 SNESSolve_NGMRES() line 262 in /Users/hall/Documents/Fortran_Codes/Packages/petsc-3.9.3/src/snes/impls/ngmres/snesngmres.c<br>
[0]PETSC ERROR: #6 SNESSolve() line 4314 in /Users/hall/Documents/Fortran_Codes/Packages/petsc-3.9.3/src/snes/interface/snes.c<br>
<br>
So it appears that Petsc is trying to do a matrix multiply.  But this is a matrix free method. How do I fix this?<br>
<br>
The relevant code is attached below.  Any help with these two issues is greatly appreciated.<br>
<br>
Code follows signature.  I have omitted non-relavent portions.<br>
Kenneth Hall<br>
<br>
<br>
<br>
=============================================================<br>
<br>
<br>
   module petsc_krylov<br>
      use petsc<br>
      use data_module<br>
      use solution_module<br>
      implicit none<br>
#include <petsc/finclude/petscsys.h><br>
#include <petsc/finclude/petsc.h><br>
      private<br>
      public                              :: gmres2<br>
      SNES                                :: snes     ! nonlinear solver<br>
      PC                                  :: pc       ! preconditioner context<br>
      KSP                                 :: ksp      ! Krylov subspace method context<br>
      SNESLineSearch                      :: linesearch<br>
      Vec                                 :: x,r      ! solution, residual vectors<br>
      PetscErrorCode                      :: ierr<br>
      PetscInt                            :: n<br>
      integer                             :: neq<br>
      integer                             :: icount,level=0<br>
   contains<br>
<br>
     .<br>
     .<br>
     .<br>
<br>
      subroutine gmres2<br>
      implicit none<br>
!<br>
!.... declare local variables<br>
      PetscInt             :: its                  ! iterations for convergence<br>
      PetscMPIInt          :: rank<br>
      PetscBool            :: snes_mf,snes_mf_set,user_precond<br>
!<br>
!.... set problem size<br>
      call get_equation_count<br>
      n = neq<br>
!<br>
!.... initialliize petsc<br>
      call PetscInitialize(PETSC_NULL_CHARACTER,ierr) ; if (ierr .ne. 0) then ;  print*,'Unable to initialize PETSc' ; endif<br>
!<br>
!.... Create nonlinear solver context<br>
      call SNESCreate(PETSC_COMM_WORLD,snes,ierr); CHKERRA(ierr)<br>
      call SNESSetType(snes,SNESNGMRES,ierr); CHKERRA(ierr)<br>
      call SNESSetTolerances(snes, 0.0, 0.0, 0.0, 20, 10000,ierr); CHKERRQ(ierr)<br>
!<br>
!.... Create vectors for solution and nonlinear function<br>
      call VecCreateSeq(PETSC_COMM_SELF,n,x,ierr); CHKERRA(ierr)<br>
      call VecDuplicate(x,r,ierr); CHKERRA(ierr)<br>
!<br>
!.... Set function evaluation routine and vector<br>
      call SNESSetFunction(snes,r,FormFunction,0,ierr); CHKERRA(ierr)<br>
!<br>
!.... Customize nonlinear solver; set runtime options<br>
      snes_mf  = PETSC_FALSE;<br>
      call PetscOptionsGetBool(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,"-snes_mf",snes_mf,snes_mf_set,ierr); CHKERRA(ierr)<br>
      if (snes_mf .eqv. PETSC_TRUE) then<br>
         call SNESGetKSP(snes,ksp,ierr); CHKERRA(ierr)<br>
         call KSPGetPC(ksp,pc,ierr); CHKERRA(ierr)<br>
         call SNESGetLineSearch(snes,linesearch,ierr); CHKERRA(ierr)<br>
         call SNESLineSearchSetType(linesearch,SNESLINESEARCHBT, ierr); CHKERRA(ierr)<br>
         user_precond  = PETSC_FALSE;<br>
         call PetscOptionsHasName(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,"-user_precond",user_precond,ierr); CHKERRA(ierr)<br>
         if (user_precond .eqv. PETSC_TRUE) then<br>
            call PCSetType(pc,PCSHELL,ierr); CHKERRA(ierr)<br>
            call PCShellSetApply(pc,MatrixFreePreconditioner,ierr); CHKERRA(ierr)<br>
         else<br>
            call PCSetType(pc,PCNONE,ierr); CHKERRA(ierr)<br>
         end if<br>
      end if<br>
!<br>
!.... set snes options<br>
      call SNESSetFromOptions(snes,ierr); CHKERRA(ierr)<br>
!<br>
!.... initialize solution<br>
      call load_x_from_u(x)<br>
!<br>
!.... solve<br>
      call SNESSolve(snes,PETSC_NULL_VEC,x,ierr)<br>
!<br>
!.... Free work space.  All PETSc objects should be destroyed when they are no longer needed.<br>
      call VecDestroy(x,ierr)<br>
      call VecDestroy(r,ierr)<br>
      call SNESDestroy(snes,ierr)<br>
      call PetscFinalize(ierr)<br>
   end<br>
!<br>
!====================================================================================================<br>
!====================================================================================================<br>
!<br>
      subroutine MatrixFreePreconditioner(x1,y1)<br>
      implicit none<br>
!<br>
!.... declare passed variables<br>
      Vec            :: x1,y1     .<br>
     .<br>
     .<br>
     .<br>
<br>
      end subroutine<br>
!<br>
!====================================================================================================<br>
!====================================================================================================<br>
!<br>
      subroutine FormFunction(snes,x1,f1,dummy,ierr)<br>
      implicit none<br>
!<br>
!.... declare passed variables<br>
      SNES           :: snes<br>
      Vec            :: x1,f1,x2<br>
      PetscErrorCode :: ierr<br>
      integer        :: dummy(*)<br>
     .<br>
     .<br>
     .<br>
      end subroutine FormFunction<br>
<br>
<br>
<br>
   end module<br>
<br>
<br>
<br>
<br>
</blockquote></div><br clear="all"><div><br></div>-- <br><div dir="ltr" class="gmail_signature" data-smartmail="gmail_signature"><div dir="ltr"><div><div dir="ltr"><div><div dir="ltr"><div>What most experimenters take for granted before they begin their experiments is infinitely more interesting than any results to which their experiments lead.<br>-- Norbert Wiener</div><div><br></div><div><a href="http://www.cse.buffalo.edu/~knepley/" target="_blank">https://www.cse.buffalo.edu/~knepley/</a><br></div></div></div></div></div></div></div></div>