Actual source code: ex40f90.F90

petsc-3.6.1 2015-08-06
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:       implicit none
 13: #include <petsc/finclude/petsc.h>

 15:       SNES             snes
 16:       PetscErrorCode   ierr
 17:       DM               da
 18:       PetscInt         ten,two,one
 19:       external         FormFunctionLocal


 22:       call PetscInitialize(PETSC_NULL_CHARACTER,ierr)

 24:       ten = 10
 25:       one = 1
 26:       two = 2

 28:       call DMDACreate2d(PETSC_COMM_WORLD,                               &
 29:      &     DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,                       &
 30:      &                DMDA_STENCIL_BOX,                                 &
 31:      &                -ten,-ten,PETSC_DECIDE,PETSC_DECIDE,two,one,      &
 32:      &                PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,da,ierr)


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

 37:       call SNESCreate(PETSC_COMM_WORLD,snes,ierr)
 38:       call SNESSetDM(snes,da,ierr)

 40:       call DMDASNESSetFunctionLocal(da,INSERT_VALUES,FormFunctionLocal,    &
 41:      &                              PETSC_NULL_OBJECT,ierr)
 42:       call SNESSetFromOptions(snes,ierr)

 44: !      Solve the nonlinear system
 45: !
 46:       call SNESSolve(snes,PETSC_NULL_OBJECT,PETSC_NULL_OBJECT,ierr)

 48:       call SNESDestroy(snes,ierr)
 49:       call DMDestroy(da,ierr)
 50:       call PetscFinalize(ierr)
 51:       end


 54:       subroutine FormFunctionLocal(in,x,f,dummy,ierr)
 55:       implicit none
 56:       PetscInt i,j,k,dummy
 57:       DMDALocalInfo in(DMDA_LOCAL_INFO_SIZE)
 58:       PetscScalar x(in(DMDA_LOCAL_INFO_DOF),                            &
 59:      &              XG_RANGE,                                           &
 60:      &              YG_RANGE)
 61:       PetscScalar f(in(DMDA_LOCAL_INFO_DOF),                            &
 62:      &             X_RANGE,                                             &
 63:      &             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