PETSc for Fortran Users#

Most of the functionality of PETSc can be obtained from Fortran programs. Make sure the suffix of all your Fortran files is .F90, not .f or .f90.

Modules and Include Files#

To use PETSc with Fortran you must use both PETSc include files and modules. At the beginning of every function and module definition you need something like

#include "petsc/finclude/petscXXX.h"
         use petscXXX

You can declare PETSc object variables using either of the following:

XXX variablename
type(tXXX) variablename

For example,

#include "petsc/finclude/petscvec.h"
         use petscvec

    Vec b
    type(tVec) x

PETSc types like PetscInt and PetscReal are simply aliases for basic Fortran types and cannot be written as type(tPetscInt)

The Fortran include files for PETSc are located in the directory $PETSC_DIR/include/petsc/finclude and the module files are located in $PETSC_DIR/$PETSC_ARCH/include

Most Fortran routines have the same names as the corresponding C versions, and PETSc command line options are fully supported. The routine arguments follow the usual Fortran conventions; the user need not worry about passing pointers or values. The calling sequences for the Fortran version are in most cases identical to the C version, except for the error checking variable discussed in Error Checking and a few routines listed in Routines with Different Fortran Interfaces.

When passing floating point numbers into PETSc Fortran subroutines, always make sure you have them marked as double precision (e.g., pass in 10.d0 instead of 10.0 or declare them as PETSc variables, e.g. PetscScalar one = 1.0). Otherwise, the compiler interprets the input as a single precision number, which can cause crashes or other mysterious problems. We highly recommend using the implicit none option at the beginning of each Fortran subroutine and declare all variables.

Error Checking#

In the Fortran version, each PETSc routine has as its final argument an integer error variable. The error code is set to be nonzero if an error has been detected; otherwise, it is zero. For example, the Fortran and C variants of KSPSolve() are given, respectively, below, where ierr denotes the PetscErrorCode error variable:

call KSPSolve(ksp, b, x, ierr) ! Fortran
ierr = KSPSolve(ksp, b, x);    // C

For proper error handling one should not use the above syntax instead one should use

PetscCall(KSPSolve(ksp, b, x, ierr))   ! Fortran subroutines
PetscCallA(KSPSolve(ksp, b, x, ierr))  ! Fortran main program
PetscCall(KSPSolve(ksp, b, x))         // C

Calling Fortran Routines from C (and C Routines from Fortran)#

Different compilers have different methods of naming Fortran routines called from C (or C routines called from Fortran). Most Fortran compilers change all the capital letters in Fortran routines to all lowercase. With some compilers, the Fortran compiler appends an underscore to the end of each Fortran routine name; for example, the Fortran routine Dabsc() would be called from C with dabsc_(). Other compilers change all the letters in Fortran routine names to capitals.

PETSc provides two macros (defined in C/C++) to help write portable code that mixes C/C++ and Fortran. They are PETSC_HAVE_FORTRAN_UNDERSCORE and PETSC_HAVE_FORTRAN_CAPS , which will be defined in the file $PETSC_DIR/$PETSC_ARCH/include/petscconf.h based on the compilers conventions. The macros are used, for example, as follows:

#if defined(PETSC_HAVE_FORTRAN_CAPS)
#define dabsc_ DABSC
#elif !defined(PETSC_HAVE_FORTRAN_UNDERSCORE)
#define dabsc_ dabsc
#endif
.....
dabsc_( &n,x,y); /* call the Fortran function */

Passing Null Pointers#

Many PETSc C functions have the option of passing a NULL argument (for example, the fifth argument of MatCreateSeqAIJ()). From Fortran, users must pass PETSC_NULL_XXX to indicate a null argument (where XXX is INTEGER, DOUBLE, CHARACTER, SCALAR, VEC, MAT, etc depending on the argument type); passing a literal 0 from Fortran in this case will crash the code. For example, when no options prefix is desired in the routine PetscOptionsGetInt(), one must use the following command in Fortran:

PetscCall(PetscOptionsGetInt(PETSC_NULL_OPTIONS, PETSC_NULL_CHARACTER, PETSC_NULL_CHARACTER, '-name', N, flg, ierr))

Matrix, Vector and IS Indices#

All matrices, vectors and IS in PETSc use zero-based indexing, regardless of whether C or Fortran is being used. For example, MatSetValues() and VecSetValues() always use zero indexing. See Basic Matrix Operations for further details.

Setting Routines#

When a function pointer (declared as external in Fortran) is passed as an argument to a PETSc function, such as the test function in KSPSetConvergenceTest(), it is assumed that this function references a routine written in the same language as the PETSc interface function that was called. For instance, if KSPSetConvergenceTest() is called from C, the test function must be a C function. Likewise, if it is called from Fortran, the test function must be (a subroutine) written in Fortran.

Compiling and Linking Fortran Programs#

See Writing C/C++ or Fortran Applications.

Routines with Different Fortran Interfaces#

The following Fortran routines differ slightly from their C counterparts; see the manual pages and previous discussion in this chapter for details:

The following functions are not supported in Fortran:

Duplicating Multiple Vectors#

The Fortran interface to VecDuplicateVecs() differs slightly from the C/C++ variant. To create n vectors of the same format as an existing vector, the user must declare a vector array, v_new of size n. Then, after VecDuplicateVecs() has been called, v_new will contain (pointers to) the new PETSc vector objects. When finished with the vectors, the user should destroy them by calling VecDestroyVecs(). For example, the following code fragment duplicates v_old to form two new vectors, v_new(1) and v_new(2).

Vec          v_old, v_new(2)
PetscInt     ierr
PetscScalar  alpha
....
PetscCall(VecDuplicateVecs(v_old, 2, v_new, ierr))
alpha = 4.3
PetscCall(VecSet(v_new(1), alpha, ierr))
alpha = 6.0
PetscCall(VecSet(v_new(2), alpha, ierr))
....
PetscCall(VecDestroyVecs(2, v_new, ierr))

Routines that Return Fortran Allocatable Arrays#

Many PETSc functions that return an array of values in C in an argument (such as ISGetIndices()) return an allocatable array in Fortran. The Fortran function names for these are suffixed with F90 as indicated below. A few routines, such as VecDuplicateVecs() discussed above, do not return a Fortran allocatable array; a large enough array must be explicitly declared before use.

C-API

Fortran-API

ISGetIndices()

ISGetIndicesF90()

ISRestoreIndices()

ISRestoreIndicesF90()

ISLocalToGlobalMappingGetIndices()

ISLocalToGlobalMappingGetIndicesF90()

ISLocalToGlobalMappingRestoreIndices()

ISLocalToGlobalMappingRestoreIndicesF90()

VecGetArray()

VecGetArrayF90()

VecRestoreArray()

VecRestoreArrayF90()

VecGetArrayRead()

VecGetArrayReadF90()

VecRestoreArrayRead()

VecRestoreArrayReadF90()

VecDuplicateVecs()

VecDuplicateVecsF90()

VecDestroyVecs()

VecDestroyVecsF90()

DMDAVecGetArray()

DMDAVecGetArrayF90()

DMDAVecRestoreArray()

DMDAVecRestoreArrayF90()

DMDAVecGetArrayRead()

DMDAVecGetArrayReadF90()

DMDAVecRestoreArrayRead()

DMDAVecRestoreArrayReadF90()

DMDAVecGetArrayWrite()

DMDAVecGetArrayWriteF90()

DMDAVecRestoreArrayWrite()

DMDAVecRestoreArrayWriteF90()

MatGetRowIJ()

MatGetRowIJF90()

MatRestoreRowIJ()

MatRestoreRowIJF90()

MatSeqAIJGetArray()

MatSeqAIJGetArrayF90()

MatSeqAIJRestoreArray()

MatSeqAIJRestoreArrayF90()

MatMPIAIJGetSeqAIJ()

MatMPIAIJGetSeqAIJF90()

MatDenseGetArray()

MatDenseGetArrayF90()

MatDenseRestoreArray()

MatDenseRestoreArrayF90()

The array arguments to these Fortran functions should be declared with forms such as

PetscScalar, pointer :: x(:)
PetscInt, pointer :: idx(:)

See the manual pages for details and pointers to example programs.

Sample Fortran Programs#

Sample programs that illustrate the PETSc interface for Fortran are given below, corresponding to Vec Test ex19f, Vec Tutorial ex4f, Draw Test ex5f, and SNES Tutorial ex1f, respectively. We also refer Fortran programmers to the C examples listed throughout the manual, since PETSc usage within the two languages differs only slightly.

Listing: src/vec/vec/tests/ex19f.F90

!
!
      program main
#include <petsc/finclude/petscvec.h>
      use petscvec
      implicit none
!
!  This example demonstrates basic use of the PETSc Fortran interface
!  to vectors.
!
       PetscInt  n
       PetscErrorCode ierr
       PetscBool  flg
       PetscScalar      one,two,three,dot
       PetscReal        norm,rdot
       Vec              x,y,w
       PetscOptions     options

       n     = 20
       one   = 1.0
       two   = 2.0
       three = 3.0

       PetscCallA(PetscInitialize(ierr))
       PetscCallA(PetscOptionsCreate(options,ierr))
       PetscCallA(PetscOptionsGetInt(options,PETSC_NULL_CHARACTER,'-n',n,flg,ierr))
       PetscCallA(PetscOptionsDestroy(options,ierr))

! Create a vector, then duplicate it
       PetscCallA(VecCreate(PETSC_COMM_WORLD,x,ierr))
       PetscCallA(VecSetSizes(x,PETSC_DECIDE,n,ierr))
       PetscCallA(VecSetFromOptions(x,ierr))
       PetscCallA(VecDuplicate(x,y,ierr))
       PetscCallA(VecDuplicate(x,w,ierr))

       PetscCallA(VecSet(x,one,ierr))
       PetscCallA(VecSet(y,two,ierr))

       PetscCallA(VecDot(x,y,dot,ierr))
       rdot = PetscRealPart(dot)
       write(6,100) rdot
  100  format('Result of inner product ',f10.4)

       PetscCallA(VecScale(x,two,ierr))
       PetscCallA(VecNorm(x,NORM_2,norm,ierr))
       write(6,110) norm
  110  format('Result of scaling ',f10.4)

       PetscCallA(VecCopy(x,w,ierr))
       PetscCallA(VecNorm(w,NORM_2,norm,ierr))
       write(6,120) norm
  120  format('Result of copy ',f10.4)

       PetscCallA(VecAXPY(y,three,x,ierr))
       PetscCallA(VecNorm(y,NORM_2,norm,ierr))
       write(6,130) norm
  130  format('Result of axpy ',f10.4)

       PetscCallA(VecDestroy(x,ierr))
       PetscCallA(VecDestroy(y,ierr))
       PetscCallA(VecDestroy(w,ierr))
       PetscCallA(PetscFinalize(ierr))
       end

Listing: src/vec/vec/tutorials/ex4f.F90

!
!
!  Description:  Illustrates the use of VecSetValues() to set
!  multiple values at once; demonstrates VecGetArrayF90().
!
! -----------------------------------------------------------------------

      program main
#include <petsc/finclude/petscvec.h>
      use petscvec
      implicit none

! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
!                 Beginning of program
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

       PetscScalar xwork(6)
       PetscScalar, pointer :: xx_v(:),yy_v(:)
       PetscInt     i,n,loc(6),isix
       PetscErrorCode ierr
       Vec         x,y

       PetscCallA(PetscInitialize(ierr))
       n = 6
       isix = 6

!  Create initial vector and duplicate it

       PetscCallA(VecCreateSeq(PETSC_COMM_SELF,n,x,ierr))
       PetscCallA(VecDuplicate(x,y,ierr))

!  Fill work arrays with vector entries and locations.  Note that
!  the vector indices are 0-based in PETSc (for both Fortran and
!  C vectors)

       do 10 i=1,n
          loc(i) = i-1
          xwork(i) = 10.0*real(i)
  10   continue

!  Set vector values.  Note that we set multiple entries at once.
!  Of course, usually one would create a work array that is the
!  natural size for a particular problem (not one that is as long
!  as the full vector).

       PetscCallA(VecSetValues(x,isix,loc,xwork,INSERT_VALUES,ierr))

!  Assemble vector

       PetscCallA(VecAssemblyBegin(x,ierr))
       PetscCallA(VecAssemblyEnd(x,ierr))

!  View vector
       PetscCallA(PetscObjectSetName(x, 'initial vector:',ierr))
       PetscCallA(VecView(x,PETSC_VIEWER_STDOUT_SELF,ierr))
       PetscCallA(VecCopy(x,y,ierr))

!  Get a pointer to vector data.
!    - For default PETSc vectors, VecGetArray() returns a pointer to
!      the data array.  Otherwise, the routine is implementation dependent.
!    - You MUST call VecRestoreArray() when you no longer need access to
!      the array.
!    - Note that the Fortran interface to VecGetArray() differs from the
!      C version.  See the users manual for details.

       PetscCallA(VecGetArrayF90(x,xx_v,ierr))
       PetscCallA(VecGetArrayF90(y,yy_v,ierr))

!  Modify vector data

       do 30 i=1,n
          xx_v(i) = 100.0*real(i)
          yy_v(i) = 1000.0*real(i)
  30   continue

!  Restore vectors

       PetscCallA(VecRestoreArrayF90(x,xx_v,ierr))
       PetscCallA(VecRestoreArrayF90(y,yy_v,ierr))

!  View vectors
       PetscCallA(PetscObjectSetName(x, 'new vector 1:',ierr))
       PetscCallA(VecView(x,PETSC_VIEWER_STDOUT_SELF,ierr))

       PetscCallA(PetscObjectSetName(y, 'new vector 2:',ierr))
       PetscCallA(VecView(y,PETSC_VIEWER_STDOUT_SELF,ierr))

!  Free work space.  All PETSc objects should be destroyed when they
!  are no longer needed.

       PetscCallA(VecDestroy(x,ierr))
       PetscCallA(VecDestroy(y,ierr))
       PetscCallA(PetscFinalize(ierr))
       end

Listing: src/sys/classes/draw/tests/ex5f.F90

!
!
      program main
#include <petsc/finclude/petscsys.h>
#include <petsc/finclude/petscdraw.h>
      use petscsys
      implicit none
!
!  This example demonstrates basic use of the Fortran interface for
!  PetscDraw routines.
!
      PetscDraw         draw
      PetscDrawLG       lg
      PetscDrawAxis     axis
      PetscErrorCode    ierr
      PetscBool         flg
      integer           x,y,width,height
      PetscScalar       xd,yd
      PetscReal         ten
      PetscInt          i,n,w,h
      PetscInt          one

      n      = 15
      x      = 0
      y      = 0
      w      = 400
      h      = 300
      ten    = 10.0
      one    = 1

      PetscCallA(PetscInitialize(ierr))

!  GetInt requires a PetscInt so have to do this ugly setting
      PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-width',w, flg,ierr))
      width = int(w)
      PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-height',h,flg,ierr))
      height = int(h)
      PetscCallA(PetscOptionsGetInt(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-n',n,flg,ierr))

      PetscCallA(PetscDrawCreate(PETSC_COMM_WORLD,PETSC_NULL_CHARACTER,PETSC_NULL_CHARACTER,x,y,width,height,draw,ierr))
      PetscCallA(PetscDrawSetFromOptions(draw,ierr))

      PetscCallA(PetscDrawLGCreate(draw,one,lg,ierr))
      PetscCallA(PetscDrawLGGetAxis(lg,axis,ierr))
      PetscCallA(PetscDrawAxisSetColors(axis,PETSC_DRAW_BLACK,PETSC_DRAW_RED,PETSC_DRAW_BLUE,ierr))
      PetscCallA(PetscDrawAxisSetLabels(axis,'toplabel','xlabel','ylabel',ierr))

      do 10, i=0,n-1
        xd = real(i) - 5.0
        yd = xd*xd
        PetscCallA(PetscDrawLGAddPoint(lg,xd,yd,ierr))
 10   continue

      PetscCallA(PetscDrawLGSetUseMarkers(lg,PETSC_TRUE,ierr))
      PetscCallA(PetscDrawLGDraw(lg,ierr))

      PetscCallA(PetscSleep(ten,ierr))

      PetscCallA(PetscDrawLGDestroy(lg,ierr))
      PetscCallA(PetscDrawDestroy(draw,ierr))
      PetscCallA(PetscFinalize(ierr))
      end

Listing: src/snes/tutorials/ex1f.F90

!
!
!  Description: Uses the Newton method to solve a two-variable system.
!

      program main
#include <petsc/finclude/petsc.h>
      use petsc
      implicit none

! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
!                   Variable declarations
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
!
!  Variables:
!     snes        - nonlinear solver
!     ksp        - linear solver
!     pc          - preconditioner context
!     ksp         - Krylov subspace method context
!     x, r        - solution, residual vectors
!     J           - Jacobian matrix
!     its         - iterations for convergence
!
      SNES     snes
      PC       pc
      KSP      ksp
      Vec      x,r
      Mat      J
      SNESLineSearch linesearch
      PetscErrorCode  ierr
      PetscInt its,i2,i20
      PetscMPIInt size,rank
      PetscScalar   pfive
      PetscReal   tol
      PetscBool   setls
#if defined(PETSC_USE_LOG)
      PetscViewer viewer
#endif
      double precision threshold,oldthreshold

!  Note: Any user-defined Fortran routines (such as FormJacobian)
!  MUST be declared as external.

      external FormFunction, FormJacobian, MyLineSearch

! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
!                 Beginning of program
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

      PetscCallA(PetscInitialize(ierr))
      PetscCallA(PetscLogNestedBegin(ierr))
      threshold = 1.0
      PetscCallA(PetscLogSetThreshold(threshold,oldthreshold,ierr))
      PetscCallMPIA(MPI_Comm_size(PETSC_COMM_WORLD,size,ierr))
      PetscCallMPIA(MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr))
      PetscCheckA(size .eq. 1,PETSC_COMM_SELF,PETSC_ERR_WRONG_MPI_SIZE,'Uniprocessor example')

      i2  = 2
      i20 = 20
! - - - - - - - - - -- - - - - - - - - - - - - - - - - - - - - - - - - -
!  Create nonlinear solver context
! - - - - - - - - - -- - - - - - - - - - - - - - - - - - - - - - - - - -

      PetscCallA(SNESCreate(PETSC_COMM_WORLD,snes,ierr))

! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
!  Create matrix and vector data structures; set corresponding routines
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

!  Create vectors for solution and nonlinear function

      PetscCallA(VecCreateSeq(PETSC_COMM_SELF,i2,x,ierr))
      PetscCallA(VecDuplicate(x,r,ierr))

!  Create Jacobian matrix data structure

      PetscCallA(MatCreate(PETSC_COMM_SELF,J,ierr))
      PetscCallA(MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,i2,i2,ierr))
      PetscCallA(MatSetFromOptions(J,ierr))
      PetscCallA(MatSetUp(J,ierr))

!  Set function evaluation routine and vector

      PetscCallA(SNESSetFunction(snes,r,FormFunction,0,ierr))

!  Set Jacobian matrix data structure and Jacobian evaluation routine

      PetscCallA(SNESSetJacobian(snes,J,J,FormJacobian,0,ierr))

! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
!  Customize nonlinear solver; set runtime options
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

!  Set linear solver defaults for this problem. By extracting the
!  KSP, KSP, and PC contexts from the SNES context, we can then
!  directly call any KSP, KSP, and PC routines to set various options.

      PetscCallA(SNESGetKSP(snes,ksp,ierr))
      PetscCallA(KSPGetPC(ksp,pc,ierr))
      PetscCallA(PCSetType(pc,PCNONE,ierr))
      tol = 1.e-4
      PetscCallA(KSPSetTolerances(ksp,tol,PETSC_DEFAULT_REAL,PETSC_DEFAULT_REAL,i20,ierr))

!  Set SNES/KSP/KSP/PC runtime options, e.g.,
!      -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc>
!  These options will override those specified above as long as
!  SNESSetFromOptions() is called _after_ any other customization
!  routines.

      PetscCallA(SNESSetFromOptions(snes,ierr))

      PetscCallA(PetscOptionsHasName(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-setls',setls,ierr))

      if (setls) then
        PetscCallA(SNESGetLineSearch(snes, linesearch, ierr))
        PetscCallA(SNESLineSearchSetType(linesearch, 'shell', ierr))
        PetscCallA(SNESLineSearchShellSetUserFunc(linesearch, MyLineSearch,0,ierr))
      endif

! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
!  Evaluate initial guess; then solve nonlinear system
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

!  Note: The user should initialize the vector, x, with the initial guess
!  for the nonlinear solver prior to calling SNESSolve().  In particular,
!  to employ an initial guess of zero, the user should explicitly set
!  this vector to zero by calling VecSet().

      pfive = 0.5
      PetscCallA(VecSet(x,pfive,ierr))
      PetscCallA(SNESSolve(snes,PETSC_NULL_VEC,x,ierr))

!  View solver converged reason; we could instead use the option -snes_converged_reason
      PetscCallA(SNESConvergedReasonView(snes,PETSC_VIEWER_STDOUT_WORLD,ierr))

      PetscCallA(SNESGetIterationNumber(snes,its,ierr))
      if (rank .eq. 0) then
         write(6,100) its
      endif
  100 format('Number of SNES iterations = ',i5)

! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
!  Free work space.  All PETSc objects should be destroyed when they
!  are no longer needed.
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -

      PetscCallA(VecDestroy(x,ierr))
      PetscCallA(VecDestroy(r,ierr))
      PetscCallA(MatDestroy(J,ierr))
      PetscCallA(SNESDestroy(snes,ierr))
#if defined(PETSC_USE_LOG)
      PetscCallA(PetscViewerASCIIOpen(PETSC_COMM_WORLD,'filename.xml',viewer,ierr))
      PetscCallA(PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_XML,ierr))
      PetscCallA(PetscLogView(viewer,ierr))
      PetscCallA(PetscViewerDestroy(viewer,ierr))
#endif
      PetscCallA(PetscFinalize(ierr))
      end
!
! ------------------------------------------------------------------------
!
!  FormFunction - Evaluates nonlinear function, F(x).
!
!  Input Parameters:
!  snes - the SNES context
!  x - input vector
!  dummy - optional user-defined context (not used here)
!
!  Output Parameter:
!  f - function vector
!
      subroutine FormFunction(snes,x,f,dummy,ierr)
      use petscsnes
      implicit none

      SNES     snes
      Vec      x,f
      PetscErrorCode ierr
      integer dummy(*)

!  Declarations for use with local arrays
      PetscScalar,pointer :: lx_v(:),lf_v(:)

!  Get pointers to vector data.
!    - VecGetArrayF90() returns a pointer to the data array.
!    - You MUST call VecRestoreArrayF90() when you no longer need access to
!      the array.

      PetscCall(VecGetArrayReadF90(x,lx_v,ierr))
      PetscCall(VecGetArrayF90(f,lf_v,ierr))

!  Compute function

      lf_v(1) = lx_v(1)*lx_v(1) + lx_v(1)*lx_v(2) - 3.0
      lf_v(2) = lx_v(1)*lx_v(2) + lx_v(2)*lx_v(2) - 6.0

!  Restore vectors

      PetscCall(VecRestoreArrayReadF90(x,lx_v,ierr))
      PetscCall(VecRestoreArrayF90(f,lf_v,ierr))

      return
      end

! ---------------------------------------------------------------------
!
!  FormJacobian - Evaluates Jacobian matrix.
!
!  Input Parameters:
!  snes - the SNES context
!  x - input vector
!  dummy - optional user-defined context (not used here)
!
!  Output Parameters:
!  A - Jacobian matrix
!  B - optionally different preconditioning matrix
!
      subroutine FormJacobian(snes,X,jac,B,dummy,ierr)
      use petscsnes
      implicit none

      SNES         snes
      Vec          X
      Mat          jac,B
      PetscScalar  A(4)
      PetscErrorCode ierr
      PetscInt idx(2),i2
      integer dummy(*)

!  Declarations for use with local arrays

      PetscScalar,pointer :: lx_v(:)

!  Get pointer to vector data

      i2 = 2
      PetscCall(VecGetArrayReadF90(x,lx_v,ierr))

!  Compute Jacobian entries and insert into matrix.
!   - Since this is such a small problem, we set all entries for
!     the matrix at once.
!   - Note that MatSetValues() uses 0-based row and column numbers
!     in Fortran as well as in C (as set here in the array idx).

      idx(1) = 0
      idx(2) = 1
      A(1) = 2.0*lx_v(1) + lx_v(2)
      A(2) = lx_v(1)
      A(3) = lx_v(2)
      A(4) = lx_v(1) + 2.0*lx_v(2)
      PetscCall(MatSetValues(B,i2,idx,i2,idx,A,INSERT_VALUES,ierr))

!  Restore vector

      PetscCall(VecRestoreArrayReadF90(x,lx_v,ierr))

!  Assemble matrix

      PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY,ierr))
      PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY,ierr))
      if (B .ne. jac) then
        PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY,ierr))
        PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY,ierr))
      endif

      return
      end

      subroutine MyLineSearch(linesearch, lctx, ierr)
      use petscsnes
      implicit none

      SNESLineSearch    linesearch
      SNES              snes
      integer           lctx
      Vec               x, f,g, y, w
      PetscReal         ynorm,gnorm,xnorm
      PetscErrorCode    ierr

      PetscScalar       mone

      mone = -1.0
      PetscCall(SNESLineSearchGetSNES(linesearch, snes, ierr))
      PetscCall(SNESLineSearchGetVecs(linesearch, x, f, y, w, g, ierr))
      PetscCall(VecNorm(y,NORM_2,ynorm,ierr))
      PetscCall(VecAXPY(x,mone,y,ierr))
      PetscCall(SNESComputeFunction(snes,x,f,ierr))
      PetscCall(VecNorm(f,NORM_2,gnorm,ierr))
      PetscCall(VecNorm(x,NORM_2,xnorm,ierr))
      PetscCall(VecNorm(y,NORM_2,ynorm,ierr))
      PetscCall(SNESLineSearchSetNorms(linesearch, xnorm, gnorm, ynorm,ierr))
      return
      end