Actual source code: ex40f90.F90

  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().
 10: #include <petsc/finclude/petscsnes.h>
 11: #include <petsc/finclude/petscdmda.h>
 12: module ex40module
 13:   use petscdmda
 14:   implicit none
 15: contains
 16:   subroutine FormFunctionLocal(in, x, f, dummy, ierr)
 17:     PetscInt, intent(in) :: dummy
 18:     DMDALocalInfo in
 19:     PetscScalar, intent(in) :: x(in%DOF, in%GXS + 1:in%GXS + in%GXM, in%GYS + 1:in%GYS + in%GYM)
 20:     PetscScalar, intent(out) :: f(in%DOF, in%XS + 1:in%XS + in%XM, in%YS + 1:in%YS + in%YM)
 21:     PetscErrorCode, intent(out) :: ierr

 23:     f = x**2 - 2.0

 25:     ierr = 0
 26:   end
 27: end module ex40module

 29: program ex40f90
 30:   use petscdmda
 31:   use petscsnes
 32:   use ex40module
 33:   implicit none

 35:   SNES snes
 36:   PetscErrorCode ierr
 37:   DM da
 38:   PetscScalar, parameter :: one = 1.0
 39:   Vec x

 41:   PetscCallA(PetscInitialize(ierr))

 43:   PetscCallA(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, 10_PETSC_INT_KIND, 10_PETSC_INT_KIND, PETSC_DECIDE, PETSC_DECIDE, 2_PETSC_INT_KIND, 1_PETSC_INT_KIND, PETSC_NULL_INTEGER_ARRAY, PETSC_NULL_INTEGER_ARRAY, da, ierr))
 44:   PetscCallA(DMSetFromOptions(da, ierr))
 45:   PetscCallA(DMSetUp(da, ierr))

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

 49:   PetscCallA(SNESCreate(PETSC_COMM_WORLD, snes, ierr))
 50:   PetscCallA(SNESSetDM(snes, da, ierr))

 52:   PetscCallA(DMDASNESSetFunctionLocal(da, INSERT_VALUES, FormFunctionLocal, 0, ierr))
 53:   PetscCallA(SNESSetFromOptions(snes, ierr))

 55: !      Solve the nonlinear system
 56: !
 57:   PetscCallA(DMCreateGlobalVector(da, x, ierr))
 58:   PetscCallA(VecSet(x, one, ierr))
 59:   PetscCallA(SNESSolve(snes, PETSC_NULL_VEC, x, ierr))

 61:   PetscCallA(VecDestroy(x, ierr))
 62:   PetscCallA(SNESDestroy(snes, ierr))
 63:   PetscCallA(DMDestroy(da, ierr))
 64:   PetscCallA(PetscFinalize(ierr))
 65: end

 67: !/*TEST
 68: !
 69: !   test:
 70: !     args: -snes_monitor_short -snes_view -da_refine 1 -pc_type mg -pc_mg_type full -ksp_type fgmres -pc_mg_galerkin pmat
 71: !     requires: !single
 72: !
 73: !TEST*/