Actual source code: ex40f90.F90

petsc-3.8.4 2018-03-24
Report Typos and Errors
  1: !
  2: !  Demonstrates use of DMDASNESSetFunctionLocal() from Fortran
  3: !
  4: !    Note: the access to the entries of the local arrays below use the Fortran
  5: !   convention of starting at zero. However calls to MatSetValues()  start at 0.
  6: !   Also note that you will have to map the i,j,k coordinates to the local PETSc ordering
  7: !   before calling MatSetValuesLocal(). Often you will find that using PETSc's default
  8: !   code for computing the Jacobian works fine and you will not need to implement
  9: !   your own FormJacobianLocal().

 11:       program ex40f90
 12:  #include <petsc/finclude/petscsnes.h>
 13:  #include <petsc/finclude/petscdmda.h>
 14:       use petscsnes
 15:       use petscdmda
 16:       implicit none

 18:       SNES             snes
 19:       PetscErrorCode   ierr
 20:       DM               da
 21:       PetscInt         ten,two,one
 22:       external         FormFunctionLocal


 25:       call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
 26:       if (ierr .ne. 0) then
 27:          print*,'PetscInitialize failed'
 28:          stop
 29:       endif

 31:       ten = 10
 32:       one = 1
 33:       two = 2

 35:       call DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,ten,ten,PETSC_DECIDE,PETSC_DECIDE,two,one, &
 36:      &     PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,da,ierr);CHKERRA(ierr)
 37:       call DMSetFromOptions(da,ierr)
 38:       call DMSetUp(da,ierr)

 40: !       Create solver object and associate it with the unknowns (on the grid)

 42:       call SNESCreate(PETSC_COMM_WORLD,snes,ierr);CHKERRA(ierr)
 43:       call SNESSetDM(snes,da,ierr);CHKERRA(ierr)

 45:       call DMDASNESSetFunctionLocal(da,INSERT_VALUES,FormFunctionLocal,0,ierr);CHKERRA(ierr)
 46:       call SNESSetFromOptions(snes,ierr);CHKERRA(ierr)

 48: !      Solve the nonlinear system
 49: !
 50:       call SNESSolve(snes,PETSC_NULL_VEC,PETSC_NULL_VEC,ierr);CHKERRA(ierr)

 52:       call SNESDestroy(snes,ierr);CHKERRA(ierr)
 53:       call DMDestroy(da,ierr);CHKERRA(ierr)
 54:       call PetscFinalize(ierr)
 55:       end


 58:       subroutine FormFunctionLocal(in,x,f,dummy,ierr)
 59:       implicit none
 60:       PetscInt i,j,k,dummy
 61:       DMDALocalInfo in(DMDA_LOCAL_INFO_SIZE)
 62:       PetscScalar x(in(DMDA_LOCAL_INFO_DOF),XG_RANGE,YG_RANGE)
 63:       PetscScalar f(in(DMDA_LOCAL_INFO_DOF),X_RANGE,Y_RANGE)
 64:       PetscErrorCode ierr

 66:       do i=in(DMDA_LOCAL_INFO_XS)+1,in(DMDA_LOCAL_INFO_XS)+in(DMDA_LOCAL_INFO_XM)
 67:          do j=in(DMDA_LOCAL_INFO_YS)+1,in(DMDA_LOCAL_INFO_YS)+in(DMDA_LOCAL_INFO_YM)
 68:             do k=1,in(DMDA_LOCAL_INFO_DOF)
 69:                f(k,i,j) = x(k,i,j)*x(k,i,j) - 2.0
 70:             enddo
 71:          enddo
 72:       enddo

 74:       return
 75:       end