Actual source code: ex12f.F
petsc-3.7.1 2016-05-15
1: !
2: !
3: ! This example demonstrates basic use of the SNES Fortran interface.
4: !
5: ! Note: The program ex10.f is the same as this example, except that it
6: ! uses the Fortran .f suffix rather than the .F suffix.
7: !
8: ! In this example the application context is a Fortran integer array:
9: ! ctx(1) = da - distributed array
10: ! 2 = F - global vector where the function is stored
11: ! 3 = xl - local work vector
12: ! 4 = comm - MPI communictor
13: ! 5 = unused
14: ! 6 = N - system size
15: !
16: ! Note: Any user-defined Fortran routines (such as FormJacobian)
17: ! MUST be declared as external.
18: !
19: !
20: ! Macros to make setting/getting values into vector clearer.
21: ! The element xx(ib) is the ibth element in the vector indicated by ctx(3)
22: #define xx(ib) vxx(ixx + (ib))
23: #define ff(ib) vff(iff + (ib))
24: #define F2(ib) vF2(iF2 + (ib))
25: program main
26: implicit none
28: #include <petsc/finclude/petscsys.h>
29: #include <petsc/finclude/petscvec.h>
30: #include <petsc/finclude/petscdm.h>
31: #include <petsc/finclude/petscdmda.h>
32: #include <petsc/finclude/petscmat.h>
33: #include <petsc/finclude/petscsnes.h>
35: PetscFortranAddr ctx(6)
36: PetscMPIInt rank,size
37: PetscErrorCode ierr
38: PetscInt N,start,end,nn,i
39: PetscInt ii,its,i1,i0,i3
40: PetscBool flg
41: SNES snes
42: Mat J
43: Vec x,r,u
44: PetscScalar xp,FF,UU,h
45: character*(10) matrixname
46: external FormJacobian,FormFunction
48: call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
49: i1 = 1
50: i0 = 0
51: i3 = 3
52: N = 10
53: call PetscOptionsGetInt(PETSC_NULL_OBJECT,PETSC_NULL_CHARACTER, &
54: & '-n',N,flg,ierr)
55: h = 1.0/(N-1.0)
56: ctx(6) = N
57: ctx(4) = PETSC_COMM_WORLD
59: call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
60: call MPI_Comm_size(PETSC_COMM_WORLD,size,ierr)
62: ! Set up data structures
63: call DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,N,i1,i1, &
64: & PETSC_NULL_INTEGER,ctx(1),ierr)
66: call DMCreateGlobalVector(ctx(1),x,ierr)
67: call DMCreateLocalVector(ctx(1),ctx(3),ierr)
69: call PetscObjectSetName(x,'Approximate Solution',ierr)
70: call VecDuplicate(x,r,ierr)
71: call VecDuplicate(x,ctx(2),ierr)
72: call VecDuplicate(x,U,ierr)
73: call PetscObjectSetName(U,'Exact Solution',ierr)
75: call MatCreateAIJ(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,N, &
76: & N,i3,PETSC_NULL_INTEGER,i0,PETSC_NULL_INTEGER,J,ierr)
77: call MatSetOption(J,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE,ierr)
78: call MatGetType(J,matrixname,ierr)
80: ! Store right-hand-side of PDE and exact solution
81: call VecGetOwnershipRange(x,start,end,ierr)
82: xp = h*start
83: nn = end - start
84: ii = start
85: do 10, i=0,nn-1
86: FF = 6.0*xp + (xp+1.e-12)**6.e0
87: UU = xp*xp*xp
88: call VecSetValues(ctx(2),i1,ii,FF,INSERT_VALUES,ierr)
89: call VecSetValues(U,i1,ii,UU,INSERT_VALUES,ierr)
90: xp = xp + h
91: ii = ii + 1
92: 10 continue
93: call VecAssemblyBegin(ctx(2),ierr)
94: call VecAssemblyEnd(ctx(2),ierr)
95: call VecAssemblyBegin(U,ierr)
96: call VecAssemblyEnd(U,ierr)
98: ! Create nonlinear solver
99: call SNESCreate(PETSC_COMM_WORLD,snes,ierr)
101: ! Set various routines and options
102: call SNESSetFunction(snes,r,FormFunction,ctx,ierr)
103: call SNESSetJacobian(snes,J,J,FormJacobian,ctx,ierr)
104: call SNESSetFromOptions(snes,ierr)
106: ! Solve nonlinear system
107: call FormInitialGuess(snes,x,ierr)
108: call SNESSolve(snes,PETSC_NULL_OBJECT,x,ierr)
109: call SNESGetIterationNumber(snes,its,ierr);
111: ! Free work space. All PETSc objects should be destroyed when they
112: ! are no longer needed.
113: call VecDestroy(x,ierr)
114: call VecDestroy(ctx(3),ierr)
115: call VecDestroy(r,ierr)
116: call VecDestroy(U,ierr)
117: call VecDestroy(ctx(2),ierr)
118: call MatDestroy(J,ierr)
119: call SNESDestroy(snes,ierr)
120: call DMDestroy(ctx(1),ierr)
121: call PetscFinalize(ierr)
122: end
125: ! -------------------- Evaluate Function F(x) ---------------------
127: subroutine FormFunction(snes,x,f,ctx,ierr)
128: implicit none
129: SNES snes
130: Vec x,f
131: PetscFortranAddr ctx(*)
132: PetscMPIInt rank,size
133: PetscInt i,s,n
134: PetscErrorCode ierr
135: PetscOffset ixx,iff,iF2
136: PetscScalar h,d,vf2(1),vxx(1),vff(1)
137: #include <petsc/finclude/petscsys.h>
138: #include <petsc/finclude/petscvec.h>
139: #include <petsc/finclude/petscdm.h>
140: #include <petsc/finclude/petscdmda.h>
141: #include <petsc/finclude/petscmat.h>
142: #include <petsc/finclude/petscsnes.h>
145: call MPI_Comm_rank(ctx(4),rank,ierr)
146: call MPI_Comm_size(ctx(4),size,ierr)
147: h = 1.0/(ctx(6) - 1.0)
148: call DMGlobalToLocalBegin(ctx(1),x,INSERT_VALUES,ctx(3),ierr)
149: call DMGlobalToLocalEnd(ctx(1),x,INSERT_VALUES,ctx(3),ierr)
151: call VecGetLocalSize(ctx(3),n,ierr)
152: if (n .gt. 1000) then
153: print*, 'Local work array not big enough'
154: call MPI_Abort(PETSC_COMM_WORLD,0,ierr)
155: endif
157: !
158: ! This sets the index ixx so that vxx(ixx+1) is the first local
159: ! element in the vector indicated by ctx(3).
160: !
161: call VecGetArrayRead(ctx(3),vxx,ixx,ierr)
162: call VecGetArray(f,vff,iff,ierr)
163: call VecGetArray(ctx(2),vF2,iF2,ierr)
165: d = h*h
167: !
168: ! Note that the array vxx() was obtained from a ghosted local vector
169: ! ctx(3) while the array vff() was obtained from the non-ghosted parallel
170: ! vector F. This is why there is a need for shift variable s. Since vff()
171: ! does not have locations for the ghost variables we need to index in it
172: ! slightly different then indexing into vxx(). For example on processor
173: ! 1 (the second processor)
174: !
175: ! xx(1) xx(2) xx(3) .....
176: ! ^^^^^^^ ^^^^^ ^^^^^
177: ! ghost value 1st local value 2nd local value
178: !
179: ! ff(1) ff(2)
180: ! ^^^^^^^ ^^^^^^^
181: ! 1st local value 2nd local value
182: !
183: if (rank .eq. 0) then
184: s = 0
185: ff(1) = xx(1)
186: else
187: s = 1
188: endif
190: do 10 i=1,n-2
191: ff(i-s+1) = d*(xx(i) - 2.0*xx(i+1) &
192: & + xx(i+2)) + xx(i+1)*xx(i+1) &
193: & - F2(i-s+1)
194: 10 continue
196: if (rank .eq. size-1) then
197: ff(n-s) = xx(n) - 1.0
198: endif
200: call VecRestoreArray(f,vff,iff,ierr)
201: call VecRestoreArrayRead(ctx(3),vxx,ixx,ierr)
202: call VecRestoreArray(ctx(2),vF2,iF2,ierr)
203: return
204: end
206: ! -------------------- Form initial approximation -----------------
208: subroutine FormInitialGuess(snes,x,ierr)
209: implicit none
210: #include <petsc/finclude/petscsys.h>
211: #include <petsc/finclude/petscvec.h>
212: #include <petsc/finclude/petscsnes.h>
213: PetscErrorCode ierr
214: Vec x
215: SNES snes
216: PetscScalar five
218: five = .5
219: call VecSet(x,five,ierr)
220: return
221: end
223: ! -------------------- Evaluate Jacobian --------------------
225: subroutine FormJacobian(snes,x,jac,B,ctx,ierr)
226: implicit none
227: #include <petsc/finclude/petscsys.h>
228: #include <petsc/finclude/petscvec.h>
229: #include <petsc/finclude/petscdm.h>
230: #include <petsc/finclude/petscdmda.h>
231: #include <petsc/finclude/petscmat.h>
232: #include <petsc/finclude/petscsnes.h>
233: SNES snes
234: Vec x
235: Mat jac,B
236: PetscFortranAddr ctx(*)
237: PetscInt ii,istart,iend
238: PetscInt i,j,n,end,start,i1
239: PetscErrorCode ierr
240: PetscMPIInt rank,size
241: PetscOffset ixx
242: PetscScalar d,A,h,vxx(1)
244: i1 = 1
245: h = 1.0/(ctx(6) - 1.0)
246: d = h*h
247: call MPI_Comm_rank(ctx(4),rank,ierr)
248: call MPI_Comm_size(ctx(4),size,ierr)
250: call VecGetArrayRead(x,vxx,ixx,ierr)
251: call VecGetOwnershipRange(x,start,end,ierr)
252: n = end - start
254: if (rank .eq. 0) then
255: A = 1.0
256: call MatSetValues(jac,i1,start,i1,start,A,INSERT_VALUES,ierr)
257: istart = 1
258: else
259: istart = 0
260: endif
261: if (rank .eq. size-1) then
262: i = INT(ctx(6)-1)
263: A = 1.0
264: call MatSetValues(jac,i1,i,i1,i,A,INSERT_VALUES,ierr)
265: iend = n-1
266: else
267: iend = n
268: endif
269: do 10 i=istart,iend-1
270: ii = i + start
271: j = start + i - 1
272: call MatSetValues(jac,i1,ii,i1,j,d,INSERT_VALUES,ierr)
273: j = start + i + 1
274: call MatSetValues(jac,i1,ii,i1,j,d,INSERT_VALUES,ierr)
275: A = -2.0*d + 2.0*xx(i+1)
276: call MatSetValues(jac,i1,ii,i1,ii,A,INSERT_VALUES,ierr)
277: 10 continue
278: call VecRestoreArrayRead(x,vxx,ixx,ierr)
279: call MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY,ierr)
280: call MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY,ierr)
281: return
282: end