Actual source code: ex12f.F
petsc-3.5.4 2015-05-23
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 <finclude/petscsys.h>
29: #include <finclude/petscvec.h>
30: #include <finclude/petscdm.h>
31: #include <finclude/petscdmda.h>
32: #include <finclude/petscmat.h>
33: #include <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_CHARACTER,'-n',N,flg,ierr)
54: h = 1.d0/(N-1.d0)
55: ctx(6) = N
56: ctx(4) = PETSC_COMM_WORLD
58: call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
59: call MPI_Comm_size(PETSC_COMM_WORLD,size,ierr)
61: ! Set up data structures
62: call DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,N,i1,i1, &
63: & PETSC_NULL_INTEGER,ctx(1),ierr)
65: call DMCreateGlobalVector(ctx(1),x,ierr)
66: call DMCreateLocalVector(ctx(1),ctx(3),ierr)
68: call PetscObjectSetName(x,'Approximate Solution',ierr)
69: call VecDuplicate(x,r,ierr)
70: call VecDuplicate(x,ctx(2),ierr)
71: call VecDuplicate(x,U,ierr)
72: call PetscObjectSetName(U,'Exact Solution',ierr)
74: call MatCreateAIJ(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,N, &
75: & N,i3,PETSC_NULL_INTEGER,i0,PETSC_NULL_INTEGER,J,ierr)
76: call MatSetOption(J,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_FALSE,ierr)
77: call MatGetType(J,matrixname,ierr)
79: ! Store right-hand-side of PDE and exact solution
80: call VecGetOwnershipRange(x,start,end,ierr)
81: xp = h*start
82: nn = end - start
83: ii = start
84: do 10, i=0,nn-1
85: FF = 6.0*xp + (xp+1.e-12)**6.e0
86: UU = xp*xp*xp
87: call VecSetValues(ctx(2),i1,ii,FF,INSERT_VALUES,ierr)
88: call VecSetValues(U,i1,ii,UU,INSERT_VALUES,ierr)
89: xp = xp + h
90: ii = ii + 1
91: 10 continue
92: call VecAssemblyBegin(ctx(2),ierr)
93: call VecAssemblyEnd(ctx(2),ierr)
94: call VecAssemblyBegin(U,ierr)
95: call VecAssemblyEnd(U,ierr)
97: ! Create nonlinear solver
98: call SNESCreate(PETSC_COMM_WORLD,snes,ierr)
100: ! Set various routines and options
101: call SNESSetFunction(snes,r,FormFunction,ctx,ierr)
102: call SNESSetJacobian(snes,J,J,FormJacobian,ctx,ierr)
103: call SNESSetFromOptions(snes,ierr)
105: ! Solve nonlinear system
106: call FormInitialGuess(snes,x,ierr)
107: call SNESSolve(snes,PETSC_NULL_OBJECT,x,ierr)
108: call SNESGetIterationNumber(snes,its,ierr);
110: ! Write results if first processor
111: if (ctx(4) .eq. 0) then
112: write(6,100) its
113: endif
114: 100 format('Number of SNES iterations = ',i5)
116: ! Free work space. All PETSc objects should be destroyed when they
117: ! are no longer needed.
118: call VecDestroy(x,ierr)
119: call VecDestroy(ctx(3),ierr)
120: call VecDestroy(r,ierr)
121: call VecDestroy(U,ierr)
122: call VecDestroy(ctx(2),ierr)
123: call MatDestroy(J,ierr)
124: call SNESDestroy(snes,ierr)
125: call DMDestroy(ctx(1),ierr)
126: call PetscFinalize(ierr)
127: end
130: ! -------------------- Evaluate Function F(x) ---------------------
132: subroutine FormFunction(snes,x,f,ctx,ierr)
133: implicit none
134: SNES snes
135: Vec x,f
136: PetscFortranAddr ctx(*)
137: PetscMPIInt rank,size
138: PetscInt i,s,n
139: PetscErrorCode ierr
140: PetscOffset ixx,iff,iF2
141: PetscScalar h,d,vf2(1),vxx(1),vff(1)
142: #include <finclude/petscsys.h>
143: #include <finclude/petscvec.h>
144: #include <finclude/petscdm.h>
145: #include <finclude/petscdmda.h>
146: #include <finclude/petscmat.h>
147: #include <finclude/petscsnes.h>
150: call MPI_Comm_rank(ctx(4),rank,ierr)
151: call MPI_Comm_size(ctx(4),size,ierr)
152: h = 1.d0/(ctx(6) - 1.d0)
153: call DMGlobalToLocalBegin(ctx(1),x,INSERT_VALUES,ctx(3),ierr)
154: call DMGlobalToLocalEnd(ctx(1),x,INSERT_VALUES,ctx(3),ierr)
156: call VecGetLocalSize(ctx(3),n,ierr)
157: if (n .gt. 1000) then
158: print*, 'Local work array not big enough'
159: call MPI_Abort(PETSC_COMM_WORLD,0,ierr)
160: endif
162: !
163: ! This sets the index ixx so that vxx(ixx+1) is the first local
164: ! element in the vector indicated by ctx(3).
165: !
166: call VecGetArray(ctx(3),vxx,ixx,ierr)
167: call VecGetArray(f,vff,iff,ierr)
168: call VecGetArray(ctx(2),vF2,iF2,ierr)
170: d = h*h
172: !
173: ! Note that the array vxx() was obtained from a ghosted local vector
174: ! ctx(3) while the array vff() was obtained from the non-ghosted parallel
175: ! vector F. This is why there is a need for shift variable s. Since vff()
176: ! does not have locations for the ghost variables we need to index in it
177: ! slightly different then indexing into vxx(). For example on processor
178: ! 1 (the second processor)
179: !
180: ! xx(1) xx(2) xx(3) .....
181: ! ^^^^^^^ ^^^^^ ^^^^^
182: ! ghost value 1st local value 2nd local value
183: !
184: ! ff(1) ff(2)
185: ! ^^^^^^^ ^^^^^^^
186: ! 1st local value 2nd local value
187: !
188: if (rank .eq. 0) then
189: s = 0
190: ff(1) = xx(1)
191: else
192: s = 1
193: endif
195: do 10 i=1,n-2
196: ff(i-s+1) = d*(xx(i) - 2.d0*xx(i+1) &
197: & + xx(i+2)) + xx(i+1)*xx(i+1) &
198: & - F2(i-s+1)
199: 10 continue
201: if (rank .eq. size-1) then
202: ff(n-s) = xx(n) - 1.d0
203: endif
205: call VecRestoreArray(f,vff,iff,ierr)
206: call VecRestoreArray(ctx(3),vxx,ixx,ierr)
207: call VecRestoreArray(ctx(2),vF2,iF2,ierr)
208: return
209: end
211: ! -------------------- Form initial approximation -----------------
213: subroutine FormInitialGuess(snes,x,ierr)
214: implicit none
215: #include <finclude/petscsys.h>
216: #include <finclude/petscvec.h>
217: #include <finclude/petscsnes.h>
218: PetscErrorCode ierr
219: Vec x
220: SNES snes
221: PetscScalar five
223: five = 5.d-1
224: call VecSet(x,five,ierr)
225: return
226: end
228: ! -------------------- Evaluate Jacobian --------------------
230: subroutine FormJacobian(snes,x,jac,B,ctx,ierr)
231: implicit none
232: #include <finclude/petscsys.h>
233: #include <finclude/petscvec.h>
234: #include <finclude/petscdm.h>
235: #include <finclude/petscdmda.h>
236: #include <finclude/petscmat.h>
237: #include <finclude/petscsnes.h>
238: SNES snes
239: Vec x
240: Mat jac,B
241: PetscFortranAddr ctx(*)
242: PetscInt ii,istart,iend
243: PetscInt i,j,n,end,start,i1
244: PetscErrorCode ierr
245: PetscMPIInt rank,size
246: PetscOffset ixx
247: PetscScalar d,A,h,vxx(1)
249: i1 = 1
250: h = 1.d0/(ctx(6) - 1.d0)
251: d = h*h
252: call MPI_Comm_rank(ctx(4),rank,ierr)
253: call MPI_Comm_size(ctx(4),size,ierr)
255: call VecGetArray(x,vxx,ixx,ierr)
256: call VecGetOwnershipRange(x,start,end,ierr)
257: n = end - start
259: if (rank .eq. 0) then
260: A = 1.0
261: call MatSetValues(jac,i1,start,i1,start,A,INSERT_VALUES,ierr)
262: istart = 1
263: else
264: istart = 0
265: endif
266: if (rank .eq. size-1) then
267: i = INT(ctx(6)-1)
268: A = 1.0
269: call MatSetValues(jac,i1,i,i1,i,A,INSERT_VALUES,ierr)
270: iend = n-1
271: else
272: iend = n
273: endif
274: do 10 i=istart,iend-1
275: ii = i + start
276: j = start + i - 1
277: call MatSetValues(jac,i1,ii,i1,j,d,INSERT_VALUES,ierr)
278: j = start + i + 1
279: call MatSetValues(jac,i1,ii,i1,j,d,INSERT_VALUES,ierr)
280: A = -2.0*d + 2.0*xx(i+1)
281: call MatSetValues(jac,i1,ii,i1,ii,A,INSERT_VALUES,ierr)
282: 10 continue
283: call VecRestoreArray(x,vxx,ixx,ierr)
284: call MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY,ierr)
285: call MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY,ierr)
286: return
287: end