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