10.2. Sample Fortran77 Programs

Up: Contents Next: Additional Information Previous: Fortran90

Sample programs that illustrate the PETSc interface for Fortran are given in Figures 15 - 18 , corresponding to ${}PETSC_DIR/src/vec/examples/tests/ex19.F, ${}PETSC_DIR/src/vec/examples/tutorials/ex4f.F, ${}PETSC_DIR/src/draw/examples/tests/ex5.F, and ${}PETSC_DIR/src/snes/examples/ex1f.F, respectively. We also refer Fortran programmers to the C examples listed throughout the manual, since PETSc usage within the two languages differs only slightly.


! 
!    "$Id: ex19.F,v 1.32 1999/01/12 23:13:42 bsmith Exp $"; 
! 
#include "include/finclude/petsc.h" 
#include "include/finclude/vec.h" 
! 
!  This example demonstrates basic use of the PETSc Fortran interface 
!  to vectors. 
! 
       integer          n, ierr,flg 
       Scalar           one, two, three, dot 
       double precision norm,rdot 
       Vec              x,y,w 
 
       n     = 20 
       one   = 1.0 
       two   = 2.0 
       three = 3.0 
 
       call PetscInitialize(PETSC_NULL_CHARACTER,ierr)                   
       call OptionsGetInt(PETSC_NULL_CHARACTER,'-n',n,flg,ierr) 
 
! Create a vector, then duplicate it 
       call VecCreate(PETSC_COMM_WORLD,PETSC_DECIDE,n,x,ierr) 
       call VecSetFromOptions(x,ierr) 
       call VecDuplicate(x,y,ierr) 
       call VecDuplicate(x,w,ierr) 
 
       call VecSet(one,x,ierr) 
       call VecSet(two,y,ierr) 
 
       call VecDot(x,y,dot,ierr) 
       rdot = PetscReal(dot) 
       write(6,100) rdot 
  100  format('Result of inner product ',f10.4) 
 
       call VecScale(two,x,ierr) 
       call VecNorm(x,NORM_2,norm,ierr) 
       write(6,110) norm 
  110  format('Result of scaling ',f10.4) 
 
       call VecCopy(x,w,ierr) 
       call VecNorm(w,NORM_2,norm,ierr) 
       write(6,120) norm 
  120  format('Result of copy ',f10.4) 
 
       call VecAXPY(three,x,y,ierr) 
       call VecNorm(y,NORM_2,norm,ierr) 
       write(6,130) norm 
  130  format('Result of axpy ',f10.4) 
 
       call VecDestroy(x,ierr) 
       call VecDestroy(y,ierr) 
       call VecDestroy(w,ierr) 
       call PetscFinalize(ierr) 
       end 
 
  

Figure 15: Sample Fortran Program: Using PETSc Vectors


! 
!      "$Id: ex4f.F,v 1.21 1998/04/15 18:00:13 balay Exp $"; 
! 
!  Description:  Illustrates the use of VecSetValues() to set 
!  multiple values at once; demonstrates VecGetArray(). 
! 
!/*T 
!   Concepts: Vectors^Assembling vectors; Using vector arrays; 
!   Routines: VecCreateSeq(); VecDuplicate(); VecSetValues(); VecView(); 
!   Routines: VecCopy(); VecView(); VecGetArray(); VecRestoreArray(); 
!   Routines: VecAssemblyBegin(); VecAssemblyEnd(); VecDestroy(); 
!   Processors: 1 
!T*/ 
! ----------------------------------------------------------------------- 
 
      program ex4f 
      implicit none 
 
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
!                    Include files 
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
! 
!  The following include statements are required for Fortran programs 
!  that use PETSc vectors: 
!     petsc.h  - base PETSc routines 
!     vec.h    - vectors 
 
#include "include/finclude/petsc.h" 
#include "include/finclude/vec.h" 
 
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -  
!                   Macro definitions 
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -  
! 
!  Macros to make clearer the process of setting values in vectors and 
!  getting values from vectors. 
! 
!   - The element xx_a(ib) is element ib+1 in the vector x 
!   - Here we add 1 to the base array index to facilitate the use of 
!     conventional Fortran 1-based array indexing. 
! 
#define xx_a(ib)  xx_v(xx_i + (ib)) 
#define yy_a(ib)  yy_v(yy_i + (ib)) 
 
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
!                 Beginning of program 
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -  
 
       Scalar      xwork(6) 
       Scalar      xx_v(1), yy_v(1) 
       integer     i, n, ierr, loc(6) 
       PetscOffset xx_i, yy_i 
       Vec         x, y 
 
       call PetscInitialize(PETSC_NULL_CHARACTER,ierr) 
       n = 6 
 
!  Create initial vector and duplicate it 
 
       call VecCreateSeq(PETSC_COMM_SELF,n,x,ierr) 
       call 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*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). 
 
       call VecSetValues(x,6,loc,xwork,INSERT_VALUES,ierr) 
 
!  Assemble vector 
 
       call VecAssemblyBegin(x,ierr) 
       call VecAssemblyEnd(x,ierr) 
 
!  View vector 
 
       write(6,20) 
  20   format('initial vector:') 
       call VecView(x,VIEWER_STDOUT_SELF,ierr) 
       call 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. 
 
       call VecGetArray(x,xx_v,xx_i,ierr) 
       call VecGetArray(y,yy_v,yy_i,ierr) 
 
!  Modify vector data 
 
       do 30 i=1,n 
          xx_a(i) = 100.0*i 
          yy_a(i) = 1000.0*i 
  30   continue 
 
!  Restore vectors 
 
       call VecRestoreArray(x,xx_v,xx_i,ierr) 
       call VecRestoreArray(y,yy_v,yy_i,ierr) 
 
!  View vectors 
 
       write(6,40) 
  40   format('new vector 1:') 
       call VecView(x,VIEWER_STDOUT_SELF,ierr) 
 
       write(6,50) 
  50   format('new vector 2:') 
       call VecView(y,VIEWER_STDOUT_SELF,ierr) 
 
!  Free work space.  All PETSc objects should be destroyed when they 
!  are no longer needed. 
 
       call VecDestroy(x,ierr) 
       call VecDestroy(y,ierr) 
       call PetscFinalize(ierr) 
       end 
  

Figure 16: Sample Fortran Program: Using VecSetValues() and VecGetArray()


! 
!    "$Id: ex5.F,v 1.21 1998/04/15 18:03:01 balay Exp $"; 
! 
#include "include/finclude/petsc.h" 
#include "include/finclude/draw.h" 
! 
!  This example demonstrates basic use of the Fortran interface for 
!  Draw routines. 
! 
      Draw              draw 
      DrawLG            lg 
      DrawAxis          axis 
      integer           n,i, ierr, x, y, width, height,flg 
      Scalar            xd,yd 
 
      n      = 20 
      x      = 0 
      y      = 0 
      width  = 300 
      height = 300 
 
      call PetscInitialize(PETSC_NULL_CHARACTER,ierr) 
 
      call OptionsGetInt(PETSC_NULL_CHARACTER,'-width',width,flg,ierr)  
      call OptionsGetInt(PETSC_NULL_CHARACTER,'-height',height,flg,ierr) 
      call OptionsGetInt(PETSC_NULL_CHARACTER,'-n',n,flg,ierr) 
 
      call DrawOpenX(PETSC_COMM_SELF,PETSC_NULL_CHARACTER,              & 
     &               PETSC_NULL_CHARACTER,x,y,width,height,draw,ierr) 
 
      call DrawLGCreate(draw,1,lg,ierr) 
      call DrawLGGetAxis(lg,axis,ierr) 
      call DrawAxisSetColors(axis,DRAW_BLACK,DRAW_RED,DRAW_BLUE,ierr) 
      call DrawAxisSetLabels(axis,'toplabel','xlabel','ylabel',ierr) 
 
      do 10, i=0,n-1 
        xd = i - 5.0 
        yd = xd*xd 
        call DrawLGAddPoint(lg,xd,yd,ierr) 
 10   continue 
 
      call DrawLGIndicateDataPoints(lg,ierr) 
      call DrawLGDraw(lg,ierr) 
      call DrawFlush(draw,ierr) 
 
      call PetscSleep(10,ierr) 
 
      call DrawLGDestroy(lg,ierr) 
      call DrawDestroy(draw,ierr) 
      call PetscFinalize(ierr) 
      end 
  

Figure 17: Sample Fortran Program: Using PETSc Draw Routines


! 
! "$Id: ex1f.F,v 1.23 1998/04/23 02:11:26 balay Exp $"; 
! 
!/*T 
!  Concepts: SNES^Solving a system of nonlinear equations (basic uniprocessor example) 
!  Routines: SNESCreate(); SNESSetFunction(); SNESSetJacobian();  
!  Routines: SNESSolve(); SNESSetFromOptions(); SNESGetSLES(); 
!  Routines: SLESGetPC(); SLESGetKSP(); KSPSetTolerances(); PCSetType(); 
!  Processors: 1 
!T*/ 
! 
!  Description: Uses the Newton method to solve a two-variable system. 
! 
! ----------------------------------------------------------------------- 
 
      program main 
      implicit none 
 
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
!                    Include files 
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
! 
!  The following include statements are generally used in SNES Fortran 
!  programs: 
!     petsc.h  - base PETSc routines 
!     vec.h    - vectors 
!     mat.h    - matrices 
!     ksp.h    - Krylov subspace methods 
!     pc.h     - preconditioners 
!     sles.h   - SLES interface 
!     snes.h   - SNES interface 
!  Other include statements may be needed if using additional PETSc 
!  routines in a Fortran program, e.g., 
!     viewer.h - viewers 
!     is.h     - index sets 
! 
#include "include/finclude/petsc.h" 
#include "include/finclude/vec.h" 
#include "include/finclude/mat.h" 
#include "include/finclude/ksp.h" 
#include "include/finclude/pc.h" 
#include "include/finclude/sles.h" 
#include "include/finclude/snes.h" 
! 
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
!                   Variable declarations 
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
! 
!  Variables: 
!     snes        - nonlinear solver 
!     sles        - linear solver 
!     pc          - preconditioner context 
!     ksp         - Krylov subspace method context 
!     x, r        - solution, residual vectors 
!     J           - Jacobian matrix 
!     its         - iterations for convergence 
! 
      SNES     snes 
      SLES     sles 
      PC       pc 
      KSP      ksp 
      Vec      x, r 
      Mat      J 
      integer  ierr, its, size, rank 
      Scalar   pfive 
      double precision   tol 
 
!  Note: Any user-defined Fortran routines (such as FormJacobian) 
!  MUST be declared as external. 
 
      external FormFunction, FormJacobian 
 
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -  
!                   Macro definitions 
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -  
! 
!  Macros to make clearer the process of setting values in vectors and 
!  getting values from vectors.  These vectors are used in the routines 
!  FormFunction() and FormJacobian(). 
!   - The element lx_a(ib) is element ib in the vector x 
! 
#define lx_a(ib) lx_v(lx_i + (ib)) 
#define lf_a(ib) lf_v(lf_i + (ib)) 
! 
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
!                 Beginning of program 
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -  
 
      call PetscInitialize(PETSC_NULL_CHARACTER,ierr) 
      call MPI_Comm_size(PETSC_COMM_WORLD,size,ierr) 
      call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr) 
      if (size .ne. 1) then 
         if (rank .eq. 0) then 
            write(6,*) 'This is a uniprocessor example only!' 
         endif 
         SETERRA(1,0,' ') 
      endif 
 
! - - - - - - - - - -- - - - - - - - - - - - - - - - - - - - - - - - - -  
!  Create nonlinear solver context 
! - - - - - - - - - -- - - - - - - - - - - - - - - - - - - - - - - - - -  
 
      call SNESCreate(PETSC_COMM_WORLD,SNES_NONLINEAR_EQUATIONS,        & 
     &     snes,ierr) 
 
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
!  Create matrix and vector data structures; set corresponding routines 
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 
!  Create vectors for solution and nonlinear function 
 
      call VecCreateSeq(PETSC_COMM_SELF,2,x,ierr) 
      call VecDuplicate(x,r,ierr) 
 
!  Create Jacobian matrix data structure 
 
      call MatCreate(PETSC_COMM_SELF,2,2,J,ierr) 
 
!  Set function evaluation routine and vector 
 
      call SNESSetFunction(snes,r,FormFunction,PETSC_NULL_OBJECT,ierr) 
 
!  Set Jacobian matrix data structure and Jacobian evaluation routine 
 
      call SNESSetJacobian(snes,J,J,FormJacobian,PETSC_NULL_OBJECT,     & 
     &     ierr) 
 
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
!  Customize nonlinear solver; set runtime options 
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 
!  Set linear solver defaults for this problem. By extracting the 
!  SLES, KSP, and PC contexts from the SNES context, we can then 
!  directly call any SLES, KSP, and PC routines to set various options. 
 
      call SNESGetSLES(snes,sles,ierr) 
      call SLESGetKSP(sles,ksp,ierr) 
      call SLESGetPC(sles,pc,ierr) 
      call PCSetType(pc,PCNONE,ierr) 
      tol = 1.e-4 
      call KSPSetTolerances(ksp,tol,PETSC_DEFAULT_DOUBLE_PRECISION,     & 
     &                      PETSC_DEFAULT_DOUBLE_PRECISION,20,ierr) 
 
!  Set SNES/SLES/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. 
 
      call SNESSetFromOptions(snes,ierr) 
 
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
!  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 
      call VecSet(pfive,x,ierr) 
      call SNESSolve(snes,x,its,ierr) 
      if (rank .eq. 0) then 
         write(6,100) its 
      endif 
  100 format('Number of Newton iterations = ',i5) 
 
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
!  Free work space.  All PETSc objects should be destroyed when they 
!  are no longer needed. 
! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 
      call VecDestroy(x,ierr) 
      call VecDestroy(r,ierr) 
      call MatDestroy(J,ierr) 
      call SNESDestroy(snes,ierr) 
      call 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) 
      implicit none 
 
#include "include/finclude/petsc.h" 
#include "include/finclude/vec.h" 
#include "include/finclude/snes.h" 
 
      SNES     snes 
      Vec      x, f 
      integer  ierr, dummy(*) 
 
!  Declarations for use with local arrays 
 
      Scalar       lx_v(1), lf_v(1) 
      PetscOffset  lx_i, lf_i  
 
!  Get pointers 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 Fortran chapter of the users manual for details. 
 
      call VecGetArray(x,lx_v,lx_i,ierr) 
      call VecGetArray(f,lf_v,lf_i,ierr) 
 
!  Compute function 
 
      lf_a(1) = lx_a(1)*lx_a(1)                                         & 
     &          + lx_a(1)*lx_a(2) - 3.0 
      lf_a(2) = lx_a(1)*lx_a(2)                                         & 
     &          + lx_a(2)*lx_a(2) - 6.0 
 
!  Restore vectors 
 
      call VecRestoreArray(x,lx_v,lx_i,ierr) 
      call VecRestoreArray(f,lf_v,lf_i,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 
!  flag - flag indicating matrix structure 
! 
      subroutine FormJacobian(snes,X,jac,B,flag,dummy) 
      implicit none 
 
#include "include/finclude/petsc.h" 
#include "include/finclude/vec.h" 
#include "include/finclude/mat.h" 
#include "include/finclude/pc.h" 
#include "include/finclude/snes.h" 
 
      SNES         snes 
      Vec          X 
      Mat          jac, B 
      MatStructure flag 
      Scalar       A(4) 
      integer      ierr, idx(2), dummy(*) 
 
!  Declarations for use with local arrays 
 
      Scalar      lx_v(1) 
      PetscOffset lx_i 
 
!  Get pointer to vector data 
 
      call VecGetArray(x,lx_v,lx_i,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_a(1) + lx_a(2) 
      A(2) = lx_a(1) 
      A(3) = lx_a(2) 
      A(4) = lx_a(1) + 2.0*lx_a(2) 
      call MatSetValues(jac,2,idx,2,idx,A,INSERT_VALUES,ierr) 
      flag = SAME_NONZERO_PATTERN 
 
!  Restore vector 
 
      call VecRestoreArray(x,lx_v,lx_i,ierr) 
 
!  Assemble matrix 
 
      call MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY,ierr) 
      call MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY,ierr) 
 
      return 
      end 
 
 

Figure 18: Sample Fortran Program: Using PETSc Nonlinear Solvers


Up: Contents Next: Additional Information Previous: Fortran90