xref: /petsc/src/snes/tutorials/ex40f90.F90 (revision 749c190bad46ba447444c173d8c7a4080c70750e) !
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>
12module ex40module
13  use petscdmda
14  implicit none
15contains
16  subroutine FormFunctionLocal(in, x, f, dummy, ierr)
17    PetscInt i, j, k, dummy
18    DMDALocalInfo in
19    PetscScalar x(in%DOF, in%GXS + 1:in%GXS + in%GXM, in%GYS + 1:in%GYS + in%GYM)
20    PetscScalar f(in%DOF, in%XS + 1:in%XS + in%XM, in%YS + 1:in%YS + in%YM)
21    PetscErrorCode ierr
22
23    do i = in%XS + 1, in%XS + in%XM
24      do j = in%YS + 1, in%YS + in%YM
25        do k = 1, in%DOF
26          f(k, i, j) = x(k, i, j)*x(k, i, j) - 2.0
27        end do
28      end do
29    end do
30
31  end
32end module ex40module
33
34program ex40f90
35  use petscdmda
36  use petscsnes
37  use ex40module
38  implicit none
39
40  SNES snes
41  PetscErrorCode ierr
42  DM da
43  PetscInt ten, two, one
44  PetscScalar sone
45  Vec x
46
47  PetscCallA(PetscInitialize(ierr))
48  ten = 10
49  one = 1
50  two = 2
51  sone = 1.0
52
53  PetscCallA(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, ten, ten, PETSC_DECIDE, PETSC_DECIDE, two, one, PETSC_NULL_INTEGER_ARRAY, PETSC_NULL_INTEGER_ARRAY, da, ierr))
54  PetscCallA(DMSetFromOptions(da, ierr))
55  PetscCallA(DMSetUp(da, ierr))
56
57!       Create solver object and associate it with the unknowns (on the grid)
58
59  PetscCallA(SNESCreate(PETSC_COMM_WORLD, snes, ierr))
60  PetscCallA(SNESSetDM(snes, da, ierr))
61
62  PetscCallA(DMDASNESSetFunctionLocal(da, INSERT_VALUES, FormFunctionLocal, 0, ierr))
63  PetscCallA(SNESSetFromOptions(snes, ierr))
64
65!      Solve the nonlinear system
66!
67  PetscCallA(DMCreateGlobalVector(da, x, ierr))
68  PetscCallA(VecSet(x, sone, ierr))
69  PetscCallA(SNESSolve(snes, PETSC_NULL_VEC, x, ierr))
70
71  PetscCallA(VecDestroy(x, ierr))
72  PetscCallA(SNESDestroy(snes, ierr))
73  PetscCallA(DMDestroy(da, ierr))
74  PetscCallA(PetscFinalize(ierr))
75end
76
77!/*TEST
78!
79!   test:
80!     args: -snes_monitor_short -snes_view -da_refine 1 -pc_type mg -pc_mg_type full -ksp_type fgmres -pc_mg_galerkin pmat
81!     requires: !single
82!
83!TEST*/
84