Actual source code: ex40f90.F90

petsc-3.3-p7 2013-05-11
  1: !
  2: !  Demonstrates use of DMDASetLocalFunction() 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 <finclude/petsc.h>

 15:       SNES             snes
 16:       PetscErrorCode   ierr
 17:       DM               da
 18:       external         FormFunctionLocal


 21:       call PetscInitialize(PETSC_NULL_CHARACTER,ierr)

 23:       call DMDACreate2d(PETSC_COMM_WORLD,                               &
 24:      &     DMDA_BOUNDARY_NONE,DMDA_BOUNDARY_NONE,                       &
 25:      &                DMDA_STENCIL_BOX,                                 &
 26:      &                -10,-10,PETSC_DECIDE,PETSC_DECIDE,2,1,            &
 27:      &                PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,da,ierr)


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

 32:       call SNESCreate(PETSC_COMM_WORLD,snes,ierr)
 33:       call SNESSetDM(snes,da,ierr)

 35:       call DMDASetLocalFunction(da,FormFunctionLocal,ierr)
 36:       call SNESSetFromOptions(snes,ierr)

 38: !      Solve the nonlinear system
 39: !
 40:       call SNESSolve(snes,PETSC_NULL_OBJECT,PETSC_NULL_OBJECT,ierr)

 42:       call SNESDestroy(snes,ierr)
 43:       call DMDestroy(da,ierr)
 44:       call PetscFinalize(ierr)
 45:       end


 48:       subroutine FormFunctionLocal(in,x,f,ierr)
 49:       implicit none
 50:       PetscInt i,j,k
 51:       DMDALocalInfo in(DMDA_LOCAL_INFO_SIZE)
 52:       PetscScalar x(in(DMDA_LOCAL_INFO_DOF),                            &
 53:      &              XG_RANGE,                                           &
 54:      &              YG_RANGE)
 55:       PetscScalar f(in(DMDA_LOCAL_INFO_DOF),                            &
 56:      &             X_RANGE,                                             &
 57:      &             Y_RANGE)
 58:       PetscErrorCode ierr

 60:       do i=in(DMDA_LOCAL_INFO_XS)+1,in(DMDA_LOCAL_INFO_XS)+in(DMDA_LOCAL_INFO_XM)
 61:          do j=in(DMDA_LOCAL_INFO_YS)+1,in(DMDA_LOCAL_INFO_YS)+in(DMDA_LOCAL_INFO_YM)
 62:             do k=1,in(DMDA_LOCAL_INFO_DOF)
 63:                f(k,i,j) = x(k,i,j)*x(k,i,j) - 2.0
 64:             enddo
 65:          enddo
 66:       enddo

 68:       return
 69:       end