Actual source code: ex1f.F
petsc-3.5.4 2015-05-23
1: !
2: ! Description: This example solves a nonlinear system on 1 processor with SNES.
3: ! We solve the Bratu (SFI - solid fuel ignition) problem in a 2D rectangular
4: ! domain. The command line options include:
5: ! -par <parameter>, where <parameter> indicates the nonlinearity of the problem
6: ! problem SFI: <parameter> = Bratu parameter (0 <= par <= 6.81)
7: ! -mx <xg>, where <xg> = number of grid points in the x-direction
8: ! -my <yg>, where <yg> = number of grid points in the y-direction
9: !
10: !/*T
11: ! Concepts: SNES^sequential Bratu example
12: ! Processors: 1
13: !T*/
14: !
15: ! --------------------------------------------------------------------------
16: !
17: ! Solid Fuel Ignition (SFI) problem. This problem is modeled by
18: ! the partial differential equation
19: !
20: ! -Laplacian u - lambda*exp(u) = 0, 0 < x,y < 1,
21: !
22: ! with boundary conditions
23: !
24: ! u = 0 for x = 0, x = 1, y = 0, y = 1.
25: !
26: ! A finite difference approximation with the usual 5-point stencil
27: ! is used to discretize the boundary value problem to obtain a nonlinear
28: ! system of equations.
29: !
30: ! The parallel version of this code is snes/examples/tutorials/ex5f.F
31: !
32: ! --------------------------------------------------------------------------
34: program main
35: implicit none
37: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
38: ! Include files
39: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
40: !
41: ! The following include statements are generally used in SNES Fortran
42: ! programs:
43: ! petscsys.h - base PETSc routines
44: ! petscvec.h - vectors
45: ! petscmat.h - matrices
46: ! petscksp.h - Krylov subspace methods
47: ! petscpc.h - preconditioners
48: ! petscsnes.h - SNES interface
49: ! In addition, we need the following for use of PETSc drawing routines
50: ! petscdraw.h - drawing routines
51: ! Other include statements may be needed if using additional PETSc
52: ! routines in a Fortran program, e.g.,
53: ! petscviewer.h - viewers
54: ! petscis.h - index sets
55: !
56: #include <finclude/petscsys.h>
57: #include <finclude/petscvec.h>
58: #include <finclude/petscis.h>
59: #include <finclude/petscdraw.h>
60: #include <finclude/petscmat.h>
61: #include <finclude/petscksp.h>
62: #include <finclude/petscpc.h>
63: #include <finclude/petscsnes.h>
64: !
65: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
66: ! Variable declarations
67: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
68: !
69: ! Variables:
70: ! snes - nonlinear solver
71: ! x, r - solution, residual vectors
72: ! J - Jacobian matrix
73: ! its - iterations for convergence
74: ! matrix_free - flag - 1 indicates matrix-free version
75: ! lambda - nonlinearity parameter
76: ! draw - drawing context
77: !
78: SNES snes
79: MatColoring mc
80: Vec x,r
81: PetscDraw draw
82: Mat J
83: PetscBool matrix_free,flg,fd_coloring
84: PetscErrorCode ierr
85: PetscInt its,N, mx,my,i5
86: PetscMPIInt size,rank
87: PetscReal lambda_max,lambda_min,lambda
88: MatFDColoring fdcoloring
89: ISColoring iscoloring
91: PetscScalar lx_v(0:1)
92: PetscOffset lx_i
94: ! Store parameters in common block
96: common /params/ lambda,mx,my
98: ! Note: Any user-defined Fortran routines (such as FormJacobian)
99: ! MUST be declared as external.
101: external FormFunction,FormInitialGuess,FormJacobian
103: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
104: ! Initialize program
105: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
107: call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
108: call MPI_Comm_size(PETSC_COMM_WORLD,size,ierr)
109: call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
111: if (size .ne. 1) then
112: if (rank .eq. 0) then
113: write(6,*) 'This is a uniprocessor example only!'
114: endif
115: SETERRQ(PETSC_COMM_SELF,1,' ',ierr)
116: endif
118: ! Initialize problem parameters
119: i5 = 5
120: lambda_max = 6.81
121: lambda_min = 0.0
122: lambda = 6.0
123: mx = 4
124: my = 4
125: call PetscOptionsGetInt(PETSC_NULL_CHARACTER,'-mx',mx,flg,ierr)
126: call PetscOptionsGetInt(PETSC_NULL_CHARACTER,'-my',my,flg,ierr)
127: call PetscOptionsGetReal(PETSC_NULL_CHARACTER,'-par',lambda, &
128: & flg,ierr)
129: if (lambda .ge. lambda_max .or. lambda .le. lambda_min) then
130: if (rank .eq. 0) write(6,*) 'Lambda is out of range'
131: SETERRQ(PETSC_COMM_SELF,1,' ',ierr)
132: endif
133: N = mx*my
135: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
136: ! Create nonlinear solver context
137: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
139: call SNESCreate(PETSC_COMM_WORLD,snes,ierr)
141: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
142: ! Create vector data structures; set function evaluation routine
143: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
145: call VecCreate(PETSC_COMM_WORLD,x,ierr)
146: call VecSetSizes(x,PETSC_DECIDE,N,ierr)
147: call VecSetFromOptions(x,ierr)
148: call VecDuplicate(x,r,ierr)
150: ! Set function evaluation routine and vector. Whenever the nonlinear
151: ! solver needs to evaluate the nonlinear function, it will call this
152: ! routine.
153: ! - Note that the final routine argument is the user-defined
154: ! context that provides application-specific data for the
155: ! function evaluation routine.
157: call SNESSetFunction(snes,r,FormFunction,PETSC_NULL_OBJECT,ierr)
159: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
160: ! Create matrix data structure; set Jacobian evaluation routine
161: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
163: ! Create matrix. Here we only approximately preallocate storage space
164: ! for the Jacobian. See the users manual for a discussion of better
165: ! techniques for preallocating matrix memory.
167: call PetscOptionsHasName(PETSC_NULL_CHARACTER,'-snes_mf', &
168: & matrix_free,ierr)
169: if (.not. matrix_free) then
170: call MatCreateSeqAIJ(PETSC_COMM_WORLD,N,N,i5,PETSC_NULL_INTEGER, &
171: & J,ierr)
172: endif
174: !
175: ! This option will cause the Jacobian to be computed via finite differences
176: ! efficiently using a coloring of the columns of the matrix.
177: !
178: fd_coloring = .false.
179: call PetscOptionsHasName(PETSC_NULL_CHARACTER,'-snes_fd_coloring', &
180: & fd_coloring,ierr)
181: if (fd_coloring) then
183: !
184: ! This initializes the nonzero structure of the Jacobian. This is artificial
185: ! because clearly if we had a routine to compute the Jacobian we won't need
186: ! to use finite differences.
187: !
188: call FormJacobian(snes,x,J,J,0,ierr)
189: !
190: ! Color the matrix, i.e. determine groups of columns that share no common
191: ! rows. These columns in the Jacobian can all be computed simulataneously.
192: !
193: call MatColoringCreate(J,mc,ierr)
194: call MatColoringSetType(mc,MATCOLORINGNATURAL,ierr)
195: call MatColoringSetFromOptions(mc,ierr)
196: call MatColoringApply(mc,iscoloring,ierr)
197: call MatColoringDestroy(mc,ierr)
198: !
199: ! Create the data structure that SNESComputeJacobianDefaultColor() uses
200: ! to compute the actual Jacobians via finite differences.
201: !
202: call MatFDColoringCreate(J,iscoloring,fdcoloring,ierr)
203: call MatFDColoringSetFunction(fdcoloring,FormFunction, &
204: & PETSC_NULL_OBJECT,ierr)
205: call MatFDColoringSetFromOptions(fdcoloring,ierr)
206: call MatFDColoringSetUp(J,iscoloring,fdcoloring,ierr)
207: !
208: ! Tell SNES to use the routine SNESComputeJacobianDefaultColor()
209: ! to compute Jacobians.
210: !
211: call SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor, &
212: & fdcoloring,ierr)
213: call ISColoringDestroy(iscoloring,ierr)
215: else if (.not. matrix_free) then
217: ! Set Jacobian matrix data structure and default Jacobian evaluation
218: ! routine. Whenever the nonlinear solver needs to compute the
219: ! Jacobian matrix, it will call this routine.
220: ! - Note that the final routine argument is the user-defined
221: ! context that provides application-specific data for the
222: ! Jacobian evaluation routine.
223: ! - The user can override with:
224: ! -snes_fd : default finite differencing approximation of Jacobian
225: ! -snes_mf : matrix-free Newton-Krylov method with no preconditioning
226: ! (unless user explicitly sets preconditioner)
227: ! -snes_mf_operator : form preconditioning matrix as set by the user,
228: ! but use matrix-free approx for Jacobian-vector
229: ! products within Newton-Krylov method
230: !
231: call SNESSetJacobian(snes,J,J,FormJacobian,PETSC_NULL_OBJECT, &
232: & ierr)
233: endif
235: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
236: ! Customize nonlinear solver; set runtime options
237: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
239: ! Set runtime options (e.g., -snes_monitor -snes_rtol <rtol> -ksp_type <type>)
241: call SNESSetFromOptions(snes,ierr)
243: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
244: ! Evaluate initial guess; then solve nonlinear system.
245: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
247: ! Note: The user should initialize the vector, x, with the initial guess
248: ! for the nonlinear solver prior to calling SNESSolve(). In particular,
249: ! to employ an initial guess of zero, the user should explicitly set
250: ! this vector to zero by calling VecSet().
252: call FormInitialGuess(x,ierr)
253: call SNESSolve(snes,PETSC_NULL_OBJECT,x,ierr)
254: call SNESGetIterationNumber(snes,its,ierr);
255: if (rank .eq. 0) then
256: write(6,100) its
257: endif
258: 100 format('Number of SNES iterations = ',i1)
260: ! PetscDraw contour plot of solution
262: call PetscDrawCreate(PETSC_COMM_WORLD,PETSC_NULL_CHARACTER, &
263: & 'Solution',300,0,300,300,draw,ierr)
264: call PetscDrawSetFromOptions(draw,ierr)
266: call VecGetArray(x,lx_v,lx_i,ierr)
267: call PetscDrawTensorContour(draw,mx,my,PETSC_NULL_REAL, &
268: & PETSC_NULL_REAL,lx_v(lx_i+1),ierr)
269: call VecRestoreArray(x,lx_v,lx_i,ierr)
271: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
272: ! Free work space. All PETSc objects should be destroyed when they
273: ! are no longer needed.
274: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
276: if (.not. matrix_free) call MatDestroy(J,ierr)
277: if (fd_coloring) call MatFDColoringDestroy(fdcoloring,ierr)
279: call VecDestroy(x,ierr)
280: call VecDestroy(r,ierr)
281: call SNESDestroy(snes,ierr)
282: call PetscDrawDestroy(draw,ierr)
283: call PetscFinalize(ierr)
284: end
286: ! ---------------------------------------------------------------------
287: !
288: ! FormInitialGuess - Forms initial approximation.
289: !
290: ! Input Parameter:
291: ! X - vector
292: !
293: ! Output Parameters:
294: ! X - vector
295: ! ierr - error code
296: !
297: ! Notes:
298: ! This routine serves as a wrapper for the lower-level routine
299: ! "ApplicationInitialGuess", where the actual computations are
300: ! done using the standard Fortran style of treating the local
301: ! vector data as a multidimensional array over the local mesh.
302: ! This routine merely accesses the local vector data via
303: ! VecGetArray() and VecRestoreArray().
304: !
305: subroutine FormInitialGuess(X,ierr)
306: implicit none
308: #include <finclude/petscsys.h>
309: #include <finclude/petscvec.h>
310: #include <finclude/petscmat.h>
311: #include <finclude/petscsnes.h>
313: ! Input/output variables:
314: Vec X
315: PetscErrorCode ierr
317: ! Declarations for use with local arrays:
318: PetscScalar lx_v(0:1)
319: PetscOffset lx_i
321: 0
323: ! Get a pointer to vector data.
324: ! - For default PETSc vectors, VecGetArray() returns a pointer to
325: ! the data array. Otherwise, the routine is implementation dependent.
326: ! - You MUST call VecRestoreArray() when you no longer need access to
327: ! the array.
328: ! - Note that the Fortran interface to VecGetArray() differs from the
329: ! C version. See the users manual for details.
331: call VecGetArray(X,lx_v,lx_i,ierr)
333: ! Compute initial guess
335: call ApplicationInitialGuess(lx_v(lx_i),ierr)
337: ! Restore vector
339: call VecRestoreArray(X,lx_v,lx_i,ierr)
341: return
342: end
344: ! ---------------------------------------------------------------------
345: !
346: ! ApplicationInitialGuess - Computes initial approximation, called by
347: ! the higher level routine FormInitialGuess().
348: !
349: ! Input Parameter:
350: ! x - local vector data
351: !
352: ! Output Parameters:
353: ! f - local vector data, f(x)
354: ! ierr - error code
355: !
356: ! Notes:
357: ! This routine uses standard Fortran-style computations over a 2-dim array.
358: !
359: subroutine ApplicationInitialGuess(x,ierr)
361: implicit none
363: ! Common blocks:
364: PetscReal lambda
365: PetscInt mx,my
366: common /params/ lambda,mx,my
368: ! Input/output variables:
369: PetscScalar x(mx,my)
370: PetscErrorCode ierr
372: ! Local variables:
373: PetscInt i,j
374: PetscScalar temp1,temp,hx,hy,one
376: ! Set parameters
378: 0
379: one = 1.0
380: hx = one/(dble(mx-1))
381: hy = one/(dble(my-1))
382: temp1 = lambda/(lambda + one)
384: do 20 j=1,my
385: temp = dble(min(j-1,my-j))*hy
386: do 10 i=1,mx
387: if (i .eq. 1 .or. j .eq. 1 &
388: & .or. i .eq. mx .or. j .eq. my) then
389: x(i,j) = 0.0
390: else
391: x(i,j) = temp1 * &
392: & sqrt(min(dble(min(i-1,mx-i)*hx),dble(temp)))
393: endif
394: 10 continue
395: 20 continue
397: return
398: end
400: ! ---------------------------------------------------------------------
401: !
402: ! FormFunction - Evaluates nonlinear function, F(x).
403: !
404: ! Input Parameters:
405: ! snes - the SNES context
406: ! X - input vector
407: ! dummy - optional user-defined context, as set by SNESSetFunction()
408: ! (not used here)
409: !
410: ! Output Parameter:
411: ! F - vector with newly computed function
412: !
413: ! Notes:
414: ! This routine serves as a wrapper for the lower-level routine
415: ! "ApplicationFunction", where the actual computations are
416: ! done using the standard Fortran style of treating the local
417: ! vector data as a multidimensional array over the local mesh.
418: ! This routine merely accesses the local vector data via
419: ! VecGetArray() and VecRestoreArray().
420: !
421: subroutine FormFunction(snes,X,F,dummy,ierr)
422: implicit none
424: #include <finclude/petscsys.h>
425: #include <finclude/petscvec.h>
426: #include <finclude/petscsnes.h>
428: ! Input/output variables:
429: SNES snes
430: Vec X,F
431: PetscFortranAddr dummy
432: PetscErrorCode ierr
434: ! Common blocks:
435: PetscReal lambda
436: PetscInt mx,my
437: common /params/ lambda,mx,my
439: ! Declarations for use with local arrays:
440: PetscScalar lx_v(0:1),lf_v(0:1)
441: PetscOffset lx_i,lf_i
443: ! Get pointers to vector data.
444: ! - For default PETSc vectors, VecGetArray() returns a pointer to
445: ! the data array. Otherwise, the routine is implementation dependent.
446: ! - You MUST call VecRestoreArray() when you no longer need access to
447: ! the array.
448: ! - Note that the Fortran interface to VecGetArray() differs from the
449: ! C version. See the Fortran chapter of the users manual for details.
451: call VecGetArray(X,lx_v,lx_i,ierr)
452: call VecGetArray(F,lf_v,lf_i,ierr)
454: ! Compute function
456: call ApplicationFunction(lx_v(lx_i),lf_v(lf_i),ierr)
458: ! Restore vectors
460: call VecRestoreArray(X,lx_v,lx_i,ierr)
461: call VecRestoreArray(F,lf_v,lf_i,ierr)
463: call PetscLogFlops(11.0d0*mx*my,ierr)
465: return
466: end
468: ! ---------------------------------------------------------------------
469: !
470: ! ApplicationFunction - Computes nonlinear function, called by
471: ! the higher level routine FormFunction().
472: !
473: ! Input Parameter:
474: ! x - local vector data
475: !
476: ! Output Parameters:
477: ! f - local vector data, f(x)
478: ! ierr - error code
479: !
480: ! Notes:
481: ! This routine uses standard Fortran-style computations over a 2-dim array.
482: !
483: subroutine ApplicationFunction(x,f,ierr)
485: implicit none
487: ! Common blocks:
488: PetscReal lambda
489: PetscInt mx,my
490: common /params/ lambda,mx,my
492: ! Input/output variables:
493: PetscScalar x(mx,my),f(mx,my)
494: PetscErrorCode ierr
496: ! Local variables:
497: PetscScalar two,one,hx,hy
498: PetscScalar hxdhy,hydhx,sc
499: PetscScalar u,uxx,uyy
500: PetscInt i,j
502: 0
503: one = 1.0
504: two = 2.0
505: hx = one/dble(mx-1)
506: hy = one/dble(my-1)
507: sc = hx*hy*lambda
508: hxdhy = hx/hy
509: hydhx = hy/hx
511: ! Compute function
513: do 20 j=1,my
514: do 10 i=1,mx
515: if (i .eq. 1 .or. j .eq. 1 &
516: & .or. i .eq. mx .or. j .eq. my) then
517: f(i,j) = x(i,j)
518: else
519: u = x(i,j)
520: uxx = hydhx * (two*u &
521: & - x(i-1,j) - x(i+1,j))
522: uyy = hxdhy * (two*u - x(i,j-1) - x(i,j+1))
523: f(i,j) = uxx + uyy - sc*exp(u)
524: endif
525: 10 continue
526: 20 continue
528: return
529: end
531: ! ---------------------------------------------------------------------
532: !
533: ! FormJacobian - Evaluates Jacobian matrix.
534: !
535: ! Input Parameters:
536: ! snes - the SNES context
537: ! x - input vector
538: ! dummy - optional user-defined context, as set by SNESSetJacobian()
539: ! (not used here)
540: !
541: ! Output Parameters:
542: ! jac - Jacobian matrix
543: ! jac_prec - optionally different preconditioning matrix (not used here)
544: ! flag - flag indicating matrix structure
545: !
546: ! Notes:
547: ! This routine serves as a wrapper for the lower-level routine
548: ! "ApplicationJacobian", where the actual computations are
549: ! done using the standard Fortran style of treating the local
550: ! vector data as a multidimensional array over the local mesh.
551: ! This routine merely accesses the local vector data via
552: ! VecGetArray() and VecRestoreArray().
553: !
554: subroutine FormJacobian(snes,X,jac,jac_prec,dummy,ierr)
555: implicit none
557: #include <finclude/petscsys.h>
558: #include <finclude/petscvec.h>
559: #include <finclude/petscmat.h>
560: #include <finclude/petscpc.h>
561: #include <finclude/petscsnes.h>
563: ! Input/output variables:
564: SNES snes
565: Vec X
566: Mat jac,jac_prec
567: PetscErrorCode ierr
568: integer dummy
570: ! Common blocks:
571: PetscReal lambda
572: PetscInt mx,my
573: common /params/ lambda,mx,my
575: ! Declarations for use with local array:
576: PetscScalar lx_v(0:1)
577: PetscOffset lx_i
579: ! Get a pointer to vector data
581: call VecGetArray(X,lx_v,lx_i,ierr)
583: ! Compute Jacobian entries
585: call ApplicationJacobian(lx_v(lx_i),jac,jac_prec,ierr)
587: ! Restore vector
589: call VecRestoreArray(X,lx_v,lx_i,ierr)
591: ! Assemble matrix
593: call MatAssemblyBegin(jac_prec,MAT_FINAL_ASSEMBLY,ierr)
594: call MatAssemblyEnd(jac_prec,MAT_FINAL_ASSEMBLY,ierr)
597: return
598: end
600: ! ---------------------------------------------------------------------
601: !
602: ! ApplicationJacobian - Computes Jacobian matrix, called by
603: ! the higher level routine FormJacobian().
604: !
605: ! Input Parameters:
606: ! x - local vector data
607: !
608: ! Output Parameters:
609: ! jac - Jacobian matrix
610: ! jac_prec - optionally different preconditioning matrix (not used here)
611: ! ierr - error code
612: !
613: ! Notes:
614: ! This routine uses standard Fortran-style computations over a 2-dim array.
615: !
616: subroutine ApplicationJacobian(x,jac,jac_prec,ierr)
617: implicit none
619: #include <finclude/petscsys.h>
620: #include <finclude/petscvec.h>
621: #include <finclude/petscmat.h>
622: #include <finclude/petscpc.h>
623: #include <finclude/petscsnes.h>
625: ! Common blocks:
626: PetscReal lambda
627: PetscInt mx,my
628: common /params/ lambda,mx,my
630: ! Input/output variables:
631: PetscScalar x(mx,my)
632: Mat jac,jac_prec
633: PetscErrorCode ierr
635: ! Local variables:
636: PetscInt i,j,row(1),col(5),i1,i5
637: PetscScalar two,one, hx,hy
638: PetscScalar hxdhy,hydhx,sc,v(5)
640: ! Set parameters
642: i1 = 1
643: i5 = 5
644: one = 1.0
645: two = 2.0
646: hx = one/dble(mx-1)
647: hy = one/dble(my-1)
648: sc = hx*hy
649: hxdhy = hx/hy
650: hydhx = hy/hx
652: ! Compute entries of the Jacobian matrix
653: ! - Here, we set all entries for a particular row at once.
654: ! - Note that MatSetValues() uses 0-based row and column numbers
655: ! in Fortran as well as in C.
657: do 20 j=1,my
658: row(1) = (j-1)*mx - 1
659: do 10 i=1,mx
660: row(1) = row(1) + 1
661: ! boundary points
662: if (i .eq. 1 .or. j .eq. 1 &
663: & .or. i .eq. mx .or. j .eq. my) then
664: call MatSetValues(jac_prec,i1,row,i1,row,one, &
665: & INSERT_VALUES,ierr)
666: ! interior grid points
667: else
668: v(1) = -hxdhy
669: v(2) = -hydhx
670: v(3) = two*(hydhx + hxdhy) &
671: & - sc*lambda*exp(x(i,j))
672: v(4) = -hydhx
673: v(5) = -hxdhy
674: col(1) = row(1) - mx
675: col(2) = row(1) - 1
676: col(3) = row(1)
677: col(4) = row(1) + 1
678: col(5) = row(1) + mx
679: call MatSetValues(jac_prec,i1,row,i5,col,v, &
680: & INSERT_VALUES,ierr)
681: endif
682: 10 continue
683: 20 continue
685: return
686: end