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().
11: program ex40f90
13: #include <petsc/finclude/petscsnes.h>
14: #include <petsc/finclude/petscdmda.h>
15: use petscsnes
16: use petscdmda
17: implicit none
19: SNES snes
20: PetscErrorCode ierr
21: DM da
22: PetscInt ten,two,one
23: external FormFunctionLocal
25: call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
26: if (ierr .ne. 0) then
27: print*,'PetscInitialize failed'
28: stop
29: endif
31: ten = 10
32: one = 1
33: two = 2
35: call DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,ten,ten,PETSC_DECIDE,PETSC_DECIDE,two,one, &
36: & PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,da,ierr);CHKERRA(ierr)
37: call DMSetFromOptions(da,ierr)
38: call DMSetUp(da,ierr)
40: ! Create solver object and associate it with the unknowns (on the grid)
42: call SNESCreate(PETSC_COMM_WORLD,snes,ierr);CHKERRA(ierr)
43: call SNESSetDM(snes,da,ierr);CHKERRA(ierr)
45: call DMDASNESSetFunctionLocal(da,INSERT_VALUES,FormFunctionLocal,0,ierr);CHKERRA(ierr)
46: call SNESSetFromOptions(snes,ierr);CHKERRA(ierr)
48: ! Solve the nonlinear system
49: !
50: call SNESSolve(snes,PETSC_NULL_VEC,PETSC_NULL_VEC,ierr);CHKERRA(ierr)
52: call SNESDestroy(snes,ierr);CHKERRA(ierr)
53: call DMDestroy(da,ierr);CHKERRA(ierr)
54: call PetscFinalize(ierr)
55: end
57: subroutine FormFunctionLocal(in,x,f,dummy,ierr)
58: implicit none
59: PetscInt i,j,k,dummy
60: DMDALocalInfo in(DMDA_LOCAL_INFO_SIZE)
61: PetscScalar x(in(DMDA_LOCAL_INFO_DOF),XG_RANGE,YG_RANGE)
62: PetscScalar f(in(DMDA_LOCAL_INFO_DOF),X_RANGE,Y_RANGE)
63: PetscErrorCode ierr
65: do i=in(DMDA_LOCAL_INFO_XS)+1,in(DMDA_LOCAL_INFO_XS)+in(DMDA_LOCAL_INFO_XM)
66: do j=in(DMDA_LOCAL_INFO_YS)+1,in(DMDA_LOCAL_INFO_YS)+in(DMDA_LOCAL_INFO_YM)
67: do k=1,in(DMDA_LOCAL_INFO_DOF)
68: f(k,i,j) = x(k,i,j)*x(k,i,j) - 2.0
69: enddo
70: enddo
71: enddo
73: return
74: end
76: !/*TEST
77: !
78: ! test:
79: ! args: -snes_monitor_short -snes_view -da_refine 1 -pc_type mg -pc_mg_type full -ksp_type fgmres -pc_mg_galerkin pmat
80: ! requires: !single
81: !
82: !TEST*/