Actual source code: ex14f.F
petsc-3.6.4 2016-04-12
1: !
2: !
3: ! This example demonstrates use of the SNES Fortran interface.
4: !
5: ! Note: The program is modified from ex12f.F
6: ! Used for testing MUMPS interface, and control when the Jacobian is rebuilt
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 = rank - processor rank
13: ! 5 = size - number of processors
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)
23: #define xx(ib) vxx(ixx + (ib))
24: #define ff(ib) vff(iff + (ib))
25: #define F2(ib) vF2(iF2 + (ib))
27: module Petscmod
28: implicit none
29: #include <petsc/finclude/petscsys.h>
30: #include <petsc/finclude/petscvec.h>
31: #include <petsc/finclude/petscvec.h90>
32: #include <petsc/finclude/petscmat.h>
33: #include <petsc/finclude/petscmat.h90>
34: #include <petsc/finclude/petscviewer.h>
35: #include <petsc/finclude/petscksp.h>
36: #include <petsc/finclude/petscpc.h>
37: #include <petsc/finclude/petscsnes.h>
38: #include <petsc/finclude/petscis.h>
39: #include <petsc/finclude/petscdm.h>
40: #include <petsc/finclude/petscdmda.h>
41: end module Petscmod
43: module Snesmonitormod
44: use Petscmod
45: implicit none
46: save
47: type monctx
48: PetscInt :: its,lag
49: end type monctx
50: end module Snesmonitormod
52: ! ---------------------------------------------------------------------
53: ! ---------------------------------------------------------------------
54: ! Subroutine FormMonitor
55: ! This function lets up keep track of the SNES progress at each step
56: ! In this routine, we determine when the Jacobian is rebuilt with the parameter 'jag'
57: !
58: ! Input Parameters:
59: ! snes - SNES nonlinear solver context
60: ! its - current nonlinear iteration, starting from a call of SNESSolve()
61: ! norm - 2-norm of current residual (may be approximate)
62: ! snesmonitor - monctx designed module (included in Snesmonitormod)
63: ! ---------------------------------------------------------------------
64: subroutine FormMonitor(snes,its,norm,snesmonitor,ierr)
66: use Snesmonitormod
67: implicit none
69: SNES :: snes
70: PetscInt :: its
71: PetscScalar :: norm
72: type(monctx) :: snesmonitor
73: PetscErrorCode :: ierr
75: c write(6,*) " "
76: c write(6,*) " its ",its,snesmonitor%its,"lag",
77: c & snesmonitor%lag
78: c call flush(6)
79: if (mod(snesmonitor%its,snesmonitor%lag).eq.0)
80: & then
81: call SNESSetLagJacobian(snes,1,ierr) ! build jacobian
82: else
83: call SNESSetLagJacobian(snes,-1,ierr) ! do NOT build jacobian
84: endif
85: snesmonitor%its = snesmonitor%its + 1
86: end subroutine FormMonitor
88: ! ---------------------------------------------------------------------
89: ! ---------------------------------------------------------------------
90: program main
92: use Petscmod
93: use Snesmonitormod
95: implicit none
97: type(monctx) :: snesmonitor
98: PetscFortranAddr ctx(6)
99: PetscMPIInt rank,size
100: PetscErrorCode ierr
101: PetscInt N,start,end,nn,i
102: PetscInt ii,its,i1,i0,i3
103: PetscBool flg
104: SNES snes
105: PC pc
106: KSP ksp
107: Mat J,F
108: Vec x,r,u
109: PetscScalar xp,FF,UU,h
110: character*(10) matrixname
111: external FormJacobian,FormFunction
112: external FormMonitor
114: call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
115: i1 = 1
116: i0 = 0
117: i3 = 3
118: N = 10
119: call PetscOptionsGetInt(PETSC_NULL_CHARACTER,'-n',N,flg,ierr)
120: h = 1.d0/(N-1.d0)
121: ctx(6) = N
123: call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
124: call MPI_Comm_size(PETSC_COMM_WORLD,size,ierr)
125: ctx(4) = rank
126: ctx(5) = size
128: ! Set up data structures
129: call DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,N,i1,i1, &
130: & PETSC_NULL_INTEGER,ctx(1),ierr)
132: call DMCreateGlobalVector(ctx(1),x,ierr)
133: call DMCreateLocalVector(ctx(1),ctx(3),ierr)
135: call PetscObjectSetName(x,'Approximate Solution',ierr)
136: call VecDuplicate(x,r,ierr)
137: call VecDuplicate(x,ctx(2),ierr)
138: call VecDuplicate(x,U,ierr)
139: call PetscObjectSetName(U,'Exact Solution',ierr)
141: call MatCreateAIJ(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,N, &
142: & N,i3,PETSC_NULL_INTEGER,i0,PETSC_NULL_INTEGER,J,ierr)
144: call MatGetType(J,matrixname,ierr)
146: ! Store right-hand-side of PDE and exact solution
147: call VecGetOwnershipRange(x,start,end,ierr)
148: xp = h*start
149: nn = end - start
150: ii = start
151: do 10, i=0,nn-1
152: FF = 6.0*xp + (xp+1.e-12)**6.e0
153: UU = xp*xp*xp
154: call VecSetValues(ctx(2),i1,ii,FF,INSERT_VALUES,ierr)
155: call VecSetValues(U,i1,ii,UU,INSERT_VALUES,ierr)
156: xp = xp + h
157: ii = ii + 1
158: 10 continue
159: call VecAssemblyBegin(ctx(2),ierr)
160: call VecAssemblyEnd(ctx(2),ierr)
161: call VecAssemblyBegin(U,ierr)
162: call VecAssemblyEnd(U,ierr)
164: ! Create nonlinear solver
165: call SNESCreate(PETSC_COMM_WORLD,snes,ierr)
167: ! Set various routines and options
168: call SNESSetFunction(snes,r,FormFunction,ctx,ierr)
169: call SNESSetJacobian(snes,J,J,FormJacobian,ctx,ierr)
170: call SNESSetLagJacobian(snes,3,ierr)
172: ! Set linear solver defaults for this problem.
173: call SNESGetKSP(snes,ksp,ierr)
174: call KSPGetPC(ksp,pc,ierr)
175: #ifdef PETSC_HAVE_MUMPS
176: call PCSetType(pc,PCLU,ierr)
177: call PCFactorSetMatSolverPackage(pc,MATSOLVERMUMPS,ierr)
178: #endif
180: snesmonitor%its = 0
181: call SNESGetLagJacobian(snes,snesmonitor%lag,ierr)
182: call SNESMonitorSet(snes,FormMonitor,snesmonitor,
183: & PETSC_NULL_FUNCTION,ierr)
184: call SNESSetFromOptions(snes,ierr)
185: call FormInitialGuess(snes,x,ierr)
187: ! Setup nonlinear solver, and call SNESSolve() for first few iterations, which calls SNESSetUp() :-(
188: !--------------------------------------------------------------------------------------------------
189: call SNESSetTolerances(snes,PETSC_DEFAULT_REAL, &
190: & PETSC_DEFAULT_REAL, &
191: & PETSC_DEFAULT_REAL, &
192: & 3,PETSC_DEFAULT_INTEGER,ierr)
193: call SNESSolve(snes,PETSC_NULL_OBJECT,x,ierr)
195: #ifdef PETSC_HAVE_MUMPS
196: ! Get PCFactor to set MUMPS matrix ordering option
197: call PCFactorGetMatrix(pc,F,ierr)
198: call MatMumpsSetIcntl(F,7,2,ierr) ! must be called before next SNESSetUp? -- not being used!!!
199: #endif
201: ! Solve nonlinear system
202: call SNESSetTolerances(snes,PETSC_DEFAULT_REAL, &
203: & PETSC_DEFAULT_REAL, &
204: & PETSC_DEFAULT_REAL, &
205: & 1000,PETSC_DEFAULT_INTEGER,ierr)
207: ! Call SNESSolve() for next few iterations
208: !--------------------------------------------------
209: snesmonitor%its = snesmonitor%its - 1 !do not count the 0-th iteration
210: call SNESSolve(snes,PETSC_NULL_OBJECT,x,ierr)
211: call SNESGetIterationNumber(snes,its,ierr);
213: ! Write results if first processor
214: if (ctx(4) .eq. 0) then
215: write(6,100) its
216: endif
217: 100 format('Number of SNES iterations = ',i5)
219: ! Free work space. All PETSc objects should be destroyed when they
220: ! are no longer needed.
221: call VecDestroy(x,ierr)
222: call VecDestroy(ctx(3),ierr)
223: call VecDestroy(r,ierr)
224: call VecDestroy(U,ierr)
225: call VecDestroy(ctx(2),ierr)
226: call MatDestroy(J,ierr)
227: call SNESDestroy(snes,ierr)
228: call DMDestroy(ctx(1),ierr)
229: call PetscFinalize(ierr)
230: end
233: ! -------------------- Evaluate Function F(x) ---------------------
235: subroutine FormFunction(snes,x,f,ctx,ierr)
236: implicit none
237: SNES snes
238: Vec x,f
239: PetscFortranAddr ctx(*)
240: PetscMPIInt rank,size
241: PetscInt i,s,n
242: PetscErrorCode ierr
243: PetscOffset ixx,iff,iF2
244: PetscScalar h,d,vf2(1),vxx(1),vff(1)
245: #include <petsc/finclude/petscsys.h>
246: #include <petsc/finclude/petscvec.h>
247: #include <petsc/finclude/petscdm.h>
248: #include <petsc/finclude/petscdmda.h>
249: #include <petsc/finclude/petscmat.h>
250: #include <petsc/finclude/petscsnes.h>
253: rank = ctx(4)
254: size = ctx(5)
255: h = 1.d0/(ctx(6) - 1.d0)
256: call DMGlobalToLocalBegin(ctx(1),x,INSERT_VALUES,ctx(3),ierr)
257: call DMGlobalToLocalEnd(ctx(1),x,INSERT_VALUES,ctx(3),ierr)
259: call VecGetLocalSize(ctx(3),n,ierr)
260: if (n .gt. 1000) then
261: print*, 'Local work array not big enough'
262: call MPI_Abort(PETSC_COMM_WORLD,0,ierr)
263: endif
265: !
266: ! This sets the index ixx so that vxx(ixx+1) is the first local
267: ! element in the vector indicated by ctx(3).
268: !
269: call VecGetArray(ctx(3),vxx,ixx,ierr)
270: call VecGetArray(f,vff,iff,ierr)
271: call VecGetArray(ctx(2),vF2,iF2,ierr)
273: d = h*h
275: !
276: ! Note that the array vxx() was obtained from a ghosted local vector
277: ! ctx(3) while the array vff() was obtained from the non-ghosted parallel
278: ! vector F. This is why there is a need for shift variable s. Since vff()
279: ! does not have locations for the ghost variables we need to index in it
280: ! slightly different then indexing into vxx(). For example on processor
281: ! 1 (the second processor)
282: !
283: ! xx(1) xx(2) xx(3) .....
284: ! ^^^^^^^ ^^^^^ ^^^^^
285: ! ghost value 1st local value 2nd local value
286: !
287: ! ff(1) ff(2)
288: ! ^^^^^^^ ^^^^^^^
289: ! 1st local value 2nd local value
290: !
291: if (rank .eq. 0) then
292: s = 0
293: ff(1) = xx(1)
294: else
295: s = 1
296: endif
298: do 10 i=1,n-2
299: ff(i-s+1) = d*(xx(i) - 2.d0*xx(i+1) &
300: & + xx(i+2)) + xx(i+1)*xx(i+1) &
301: & - F2(i-s+1)
302: 10 continue
304: if (rank .eq. size-1) then
305: ff(n-s) = xx(n) - 1.d0
306: endif
308: call VecRestoreArray(f,vff,iff,ierr)
309: call VecRestoreArray(ctx(3),vxx,ixx,ierr)
310: call VecRestoreArray(ctx(2),vF2,iF2,ierr)
311: return
312: end
314: ! -------------------- Form initial approximation -----------------
316: subroutine FormInitialGuess(snes,x,ierr)
317: implicit none
318: #include <petsc/finclude/petscsys.h>
319: #include <petsc/finclude/petscvec.h>
320: #include <petsc/finclude/petscsnes.h>
321: PetscErrorCode ierr
322: Vec x
323: SNES snes
324: PetscScalar five
326: five = 5.d-1
327: call VecSet(x,five,ierr)
328: return
329: end
331: ! -------------------- Evaluate Jacobian --------------------
333: subroutine FormJacobian(snes,x,jac,B,flag,ctx,ierr)
335: use Petscmod
336: implicit none
338: SNES snes
339: Vec x
340: Mat jac,B
341: PetscFortranAddr ctx(*)
342: PetscBool flag
343: PetscInt ii,istart,iend
344: PetscInt i,j,n,end,start,i1
345: PetscErrorCode ierr
346: PetscMPIInt rank,size
347: PetscOffset ixx
348: PetscScalar d,A,h,vxx(1)
350: rank = ctx(4)
351: size = ctx(5)
352: if (rank .eq. 0) then
353: write(6,*)" Jacobian is built ..."
354: call flush(6)
355: endif
356: i1 = 1
357: h = 1.d0/(ctx(6) - 1.d0)
358: d = h*h
359: rank = ctx(4)
360: size = ctx(5)
362: call VecGetArray(x,vxx,ixx,ierr)
363: call VecGetOwnershipRange(x,start,end,ierr)
364: n = end - start
366: if (rank .eq. 0) then
367: A = 1.0
368: call MatSetValues(jac,i1,start,i1,start,A,INSERT_VALUES,ierr)
369: istart = 1
370: else
371: istart = 0
372: endif
373: if (rank .eq. size-1) then
374: i = ctx(6)-1
375: A = 1.0
376: call MatSetValues(jac,i1,i,i1,i,A,INSERT_VALUES,ierr)
377: iend = n-1
378: else
379: iend = n
380: endif
381: do 10 i=istart,iend-1
382: ii = i + start
383: j = start + i - 1
384: call MatSetValues(jac,i1,ii,i1,j,d,INSERT_VALUES,ierr)
385: j = start + i + 1
386: call MatSetValues(jac,i1,ii,i1,j,d,INSERT_VALUES,ierr)
387: A = -2.0*d + 2.0*xx(i+1)
388: call MatSetValues(jac,i1,ii,i1,ii,A,INSERT_VALUES,ierr)
389: 10 continue
390: call VecRestoreArray(x,vxx,ixx,ierr)
391: call MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY,ierr)
392: call MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY,ierr)
393: return
394: end