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