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: PetscCallA(PetscInitialize(ierr))
26: ten = 10
27: one = 1
28: two = 2
30: 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,PETSC_NULL_INTEGER,da,ierr))
31: PetscCallA(DMSetFromOptions(da,ierr))
32: PetscCallA(DMSetUp(da,ierr))
34: ! Create solver object and associate it with the unknowns (on the grid)
36: PetscCallA(SNESCreate(PETSC_COMM_WORLD,snes,ierr))
37: PetscCallA(SNESSetDM(snes,da,ierr))
39: PetscCallA(DMDASNESSetFunctionLocal(da,INSERT_VALUES,FormFunctionLocal,0,ierr))
40: PetscCallA(SNESSetFromOptions(snes,ierr))
42: ! Solve the nonlinear system
43: !
44: PetscCallA(SNESSolve(snes,PETSC_NULL_VEC,PETSC_NULL_VEC,ierr))
46: PetscCallA(SNESDestroy(snes,ierr))
47: PetscCallA(DMDestroy(da,ierr))
48: PetscCallA(PetscFinalize(ierr))
49: end
51: subroutine FormFunctionLocal(in,x,f,dummy,ierr)
52: implicit none
53: PetscInt i,j,k,dummy
54: DMDALocalInfo in(DMDA_LOCAL_INFO_SIZE)
55: PetscScalar x(in(DMDA_LOCAL_INFO_DOF),XG_RANGE,YG_RANGE)
56: PetscScalar f(in(DMDA_LOCAL_INFO_DOF),X_RANGE,Y_RANGE)
57: PetscErrorCode ierr
59: do i=in(DMDA_LOCAL_INFO_XS)+1,in(DMDA_LOCAL_INFO_XS)+in(DMDA_LOCAL_INFO_XM)
60: do j=in(DMDA_LOCAL_INFO_YS)+1,in(DMDA_LOCAL_INFO_YS)+in(DMDA_LOCAL_INFO_YM)
61: do k=1,in(DMDA_LOCAL_INFO_DOF)
62: f(k,i,j) = x(k,i,j)*x(k,i,j) - 2.0
63: enddo
64: enddo
65: enddo
67: return
68: end
70: !/*TEST
71: !
72: ! test:
73: ! args: -snes_monitor_short -snes_view -da_refine 1 -pc_type mg -pc_mg_type full -ksp_type fgmres -pc_mg_galerkin pmat
74: ! requires: !single
75: !
76: !TEST*/