Actual source code: ex1f.F90
petsc-3.8.4 2018-03-24
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: #include <petsc/finclude/petscdraw.h>
36: #include <petsc/finclude/petscsnes.h>
37: use petscsnes
38: implicit none
39: #if defined(PETSC_USING_F90) && !defined(PETSC_USE_FORTRANKIND)
40: external SNESCOMPUTEJACOBIANDEFAULTCOLOR
41: #endif
42: !
43: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
44: ! Variable declarations
45: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
46: !
47: ! Variables:
48: ! snes - nonlinear solver
49: ! x, r - solution, residual vectors
50: ! J - Jacobian matrix
51: ! its - iterations for convergence
52: ! matrix_free - flag - 1 indicates matrix-free version
53: ! lambda - nonlinearity parameter
54: ! draw - drawing context
55: !
56: SNES snes
57: MatColoring mc
58: Vec x,r
59: PetscDraw draw
60: Mat J
61: PetscBool matrix_free,flg,fd_coloring
62: PetscErrorCode ierr
63: PetscInt its,N, mx,my,i5
64: PetscMPIInt size,rank
65: PetscReal lambda_max,lambda_min,lambda
66: MatFDColoring fdcoloring
67: ISColoring iscoloring
69: PetscScalar lx_v(0:1)
70: PetscOffset lx_i
72: ! Store parameters in common block
74: common /params/ lambda,mx,my,fdcoloring,fd_coloring
76: ! Note: Any user-defined Fortran routines (such as FormJacobian)
77: ! MUST be declared as external.
79: external FormFunction,FormInitialGuess,FormJacobian
81: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
82: ! Initialize program
83: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
85: call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
86: if (ierr .ne. 0) then
87: print*,'Unable to initialize PETSc'
88: stop
89: endif
90: call MPI_Comm_size(PETSC_COMM_WORLD,size,ierr)
91: call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
93: if (size .ne. 1) then
94: if (rank .eq. 0) then
95: write(6,*) 'This is a uniprocessor example only!'
96: endif
97: SETERRA(PETSC_COMM_SELF,1,' ')
98: endif
100: ! Initialize problem parameters
101: i5 = 5
102: lambda_max = 6.81
103: lambda_min = 0.0
104: lambda = 6.0
105: mx = 4
106: my = 4
107: call PetscOptionsGetInt(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-mx',mx,flg,ierr)
108: call PetscOptionsGetInt(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-my',my,flg,ierr)
109: call PetscOptionsGetReal(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-par',lambda,flg,ierr)
110: if (lambda .ge. lambda_max .or. lambda .le. lambda_min) then
111: if (rank .eq. 0) write(6,*) 'Lambda is out of range'
112: SETERRA(PETSC_COMM_SELF,1,' ')
113: endif
114: N = mx*my
116: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
117: ! Create nonlinear solver context
118: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
120: call SNESCreate(PETSC_COMM_WORLD,snes,ierr)
122: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
123: ! Create vector data structures; set function evaluation routine
124: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
126: call VecCreate(PETSC_COMM_WORLD,x,ierr)
127: call VecSetSizes(x,PETSC_DECIDE,N,ierr)
128: call VecSetFromOptions(x,ierr)
129: call VecDuplicate(x,r,ierr)
131: ! Set function evaluation routine and vector. Whenever the nonlinear
132: ! solver needs to evaluate the nonlinear function, it will call this
133: ! routine.
134: ! - Note that the final routine argument is the user-defined
135: ! context that provides application-specific data for the
136: ! function evaluation routine.
138: call SNESSetFunction(snes,r,FormFunction,PETSC_NULL_VEC,ierr)
140: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
141: ! Create matrix data structure; set Jacobian evaluation routine
142: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
144: ! Create matrix. Here we only approximately preallocate storage space
145: ! for the Jacobian. See the users manual for a discussion of better
146: ! techniques for preallocating matrix memory.
148: call PetscOptionsHasName(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-snes_mf',matrix_free,ierr)
149: if (.not. matrix_free) then
150: call MatCreateSeqAIJ(PETSC_COMM_WORLD,N,N,i5,PETSC_NULL_INTEGER,J,ierr)
151: endif
153: !
154: ! This option will cause the Jacobian to be computed via finite differences
155: ! efficiently using a coloring of the columns of the matrix.
156: !
157: fd_coloring = .false.
158: call PetscOptionsHasName(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER,'-snes_fd_coloring',fd_coloring,ierr)
159: if (fd_coloring) then
161: !
162: ! This initializes the nonzero structure of the Jacobian. This is artificial
163: ! because clearly if we had a routine to compute the Jacobian we won't need
164: ! to use finite differences.
165: !
166: call FormJacobian(snes,x,J,J,0,ierr)
167: !
168: ! Color the matrix, i.e. determine groups of columns that share no common
169: ! rows. These columns in the Jacobian can all be computed simulataneously.
170: !
171: call MatColoringCreate(J,mc,ierr)
172: call MatColoringSetType(mc,MATCOLORINGNATURAL,ierr)
173: call MatColoringSetFromOptions(mc,ierr)
174: call MatColoringApply(mc,iscoloring,ierr)
175: call MatColoringDestroy(mc,ierr)
176: !
177: ! Create the data structure that SNESComputeJacobianDefaultColor() uses
178: ! to compute the actual Jacobians via finite differences.
179: !
180: call MatFDColoringCreate(J,iscoloring,fdcoloring,ierr)
181: call MatFDColoringSetFunction(fdcoloring,FormFunction,0,ierr)
182: call MatFDColoringSetFromOptions(fdcoloring,ierr)
183: call MatFDColoringSetUp(J,iscoloring,fdcoloring,ierr)
184: !
185: ! Tell SNES to use the routine SNESComputeJacobianDefaultColor()
186: ! to compute Jacobians.
187: !
188: call SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,fdcoloring,ierr)
189: call ISColoringDestroy(iscoloring,ierr)
191: else if (.not. matrix_free) then
193: ! Set Jacobian matrix data structure and default Jacobian evaluation
194: ! routine. Whenever the nonlinear solver needs to compute the
195: ! Jacobian matrix, it will call this routine.
196: ! - Note that the final routine argument is the user-defined
197: ! context that provides application-specific data for the
198: ! Jacobian evaluation routine.
199: ! - The user can override with:
200: ! -snes_fd : default finite differencing approximation of Jacobian
201: ! -snes_mf : matrix-free Newton-Krylov method with no preconditioning
202: ! (unless user explicitly sets preconditioner)
203: ! -snes_mf_operator : form preconditioning matrix as set by the user,
204: ! but use matrix-free approx for Jacobian-vector
205: ! products within Newton-Krylov method
206: !
207: call SNESSetJacobian(snes,J,J,FormJacobian,0,ierr)
208: endif
210: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
211: ! Customize nonlinear solver; set runtime options
212: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
214: ! Set runtime options (e.g., -snes_monitor -snes_rtol <rtol> -ksp_type <type>)
216: call SNESSetFromOptions(snes,ierr)
218: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
219: ! Evaluate initial guess; then solve nonlinear system.
220: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
222: ! Note: The user should initialize the vector, x, with the initial guess
223: ! for the nonlinear solver prior to calling SNESSolve(). In particular,
224: ! to employ an initial guess of zero, the user should explicitly set
225: ! this vector to zero by calling VecSet().
227: call FormInitialGuess(x,ierr)
228: call SNESSolve(snes,PETSC_NULL_VEC,x,ierr)
229: call SNESGetIterationNumber(snes,its,ierr);
230: if (rank .eq. 0) then
231: write(6,100) its
232: endif
233: 100 format('Number of SNES iterations = ',i1)
235: ! PetscDraw contour plot of solution
237: call PetscDrawCreate(PETSC_COMM_WORLD,PETSC_NULL_CHARACTER,'Solution',300,0,300,300,draw,ierr)
238: call PetscDrawSetFromOptions(draw,ierr)
240: call VecGetArrayRead(x,lx_v,lx_i,ierr)
241: call PetscDrawTensorContour(draw,mx,my,PETSC_NULL_REAL,PETSC_NULL_REAL,lx_v(lx_i+1),ierr)
242: call VecRestoreArrayRead(x,lx_v,lx_i,ierr)
244: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
245: ! Free work space. All PETSc objects should be destroyed when they
246: ! are no longer needed.
247: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
249: if (.not. matrix_free) call MatDestroy(J,ierr)
250: if (fd_coloring) call MatFDColoringDestroy(fdcoloring,ierr)
252: call VecDestroy(x,ierr)
253: call VecDestroy(r,ierr)
254: call SNESDestroy(snes,ierr)
255: call PetscDrawDestroy(draw,ierr)
256: call PetscFinalize(ierr)
257: end
259: ! ---------------------------------------------------------------------
260: !
261: ! FormInitialGuess - Forms initial approximation.
262: !
263: ! Input Parameter:
264: ! X - vector
265: !
266: ! Output Parameters:
267: ! X - vector
268: ! ierr - error code
269: !
270: ! Notes:
271: ! This routine serves as a wrapper for the lower-level routine
272: ! "ApplicationInitialGuess", where the actual computations are
273: ! done using the standard Fortran style of treating the local
274: ! vector data as a multidimensional array over the local mesh.
275: ! This routine merely accesses the local vector data via
276: ! VecGetArray() and VecRestoreArray().
277: !
278: subroutine FormInitialGuess(X,ierr)
279: use petscsnes
280: implicit none
282: ! Input/output variables:
283: Vec X
284: PetscErrorCode ierr
286: ! Declarations for use with local arrays:
287: PetscScalar lx_v(0:1)
288: PetscOffset lx_i
290: 0
292: ! Get a pointer to vector data.
293: ! - For default PETSc vectors, VecGetArray() returns a pointer to
294: ! the data array. Otherwise, the routine is implementation dependent.
295: ! - You MUST call VecRestoreArray() when you no longer need access to
296: ! the array.
297: ! - Note that the Fortran interface to VecGetArray() differs from the
298: ! C version. See the users manual for details.
300: call VecGetArray(X,lx_v,lx_i,ierr)
302: ! Compute initial guess
304: call ApplicationInitialGuess(lx_v(lx_i),ierr)
306: ! Restore vector
308: call VecRestoreArray(X,lx_v,lx_i,ierr)
310: return
311: end
313: ! ---------------------------------------------------------------------
314: !
315: ! ApplicationInitialGuess - Computes initial approximation, called by
316: ! the higher level routine FormInitialGuess().
317: !
318: ! Input Parameter:
319: ! x - local vector data
320: !
321: ! Output Parameters:
322: ! f - local vector data, f(x)
323: ! ierr - error code
324: !
325: ! Notes:
326: ! This routine uses standard Fortran-style computations over a 2-dim array.
327: !
328: subroutine ApplicationInitialGuess(x,ierr)
329: use petscksp
330: implicit none
332: ! Common blocks:
333: PetscReal lambda
334: PetscInt mx,my
335: MatFDColoring fdcoloring
336: PetscBool fd_coloring
337: common /params/ lambda,mx,my,fdcoloring,fd_coloring
339: ! Input/output variables:
340: PetscScalar x(mx,my)
341: PetscErrorCode ierr
343: ! Local variables:
344: PetscInt i,j
345: PetscReal temp1,temp,hx,hy,one
347: ! Set parameters
349: 0
350: one = 1.0
351: hx = one/(mx-1)
352: hy = one/(my-1)
353: temp1 = lambda/(lambda + one)
355: do 20 j=1,my
356: temp = min(j-1,my-j)*hy
357: do 10 i=1,mx
358: if (i .eq. 1 .or. j .eq. 1 .or. i .eq. mx .or. j .eq. my) then
359: x(i,j) = 0.0
360: else
361: x(i,j) = temp1 * sqrt(min(min(i-1,mx-i)*hx,temp))
362: endif
363: 10 continue
364: 20 continue
366: return
367: end
369: ! ---------------------------------------------------------------------
370: !
371: ! FormFunction - Evaluates nonlinear function, F(x).
372: !
373: ! Input Parameters:
374: ! snes - the SNES context
375: ! X - input vector
376: ! dummy - optional user-defined context, as set by SNESSetFunction()
377: ! (not used here)
378: !
379: ! Output Parameter:
380: ! F - vector with newly computed function
381: !
382: ! Notes:
383: ! This routine serves as a wrapper for the lower-level routine
384: ! "ApplicationFunction", where the actual computations are
385: ! done using the standard Fortran style of treating the local
386: ! vector data as a multidimensional array over the local mesh.
387: ! This routine merely accesses the local vector data via
388: ! VecGetArray() and VecRestoreArray().
389: !
390: subroutine FormFunction(snes,X,F,dummy,ierr)
391: use petscsnes
392: implicit none
394: ! Input/output variables:
395: SNES snes
396: Vec X,F
397: PetscFortranAddr dummy
398: PetscErrorCode ierr
400: ! Common blocks:
401: PetscReal lambda
402: PetscInt mx,my
403: MatFDColoring fdcoloring
404: PetscBool fd_coloring
405: common /params/ lambda,mx,my,fdcoloring,fd_coloring
407: ! Declarations for use with local arrays:
408: PetscScalar lx_v(0:1),lf_v(0:1)
409: PetscOffset lx_i,lf_i
411: PetscInt, pointer :: indices(:)
413: ! Get pointers to vector data.
414: ! - For default PETSc vectors, VecGetArray() returns a pointer to
415: ! the data array. Otherwise, the routine is implementation dependent.
416: ! - You MUST call VecRestoreArray() when you no longer need access to
417: ! the array.
418: ! - Note that the Fortran interface to VecGetArray() differs from the
419: ! C version. See the Fortran chapter of the users manual for details.
421: call VecGetArrayRead(X,lx_v,lx_i,ierr)
422: call VecGetArray(F,lf_v,lf_i,ierr)
424: ! Compute function
426: call ApplicationFunction(lx_v(lx_i),lf_v(lf_i),ierr)
428: ! Restore vectors
430: call VecRestoreArrayRead(X,lx_v,lx_i,ierr)
431: call VecRestoreArray(F,lf_v,lf_i,ierr)
433: call PetscLogFlops(11.0d0*mx*my,ierr)
434: !
435: ! fdcoloring is in the common block and used here ONLY to test the
436: ! calls to MatFDColoringGetPerturbedColumnsF90() and MatFDColoringRestorePerturbedColumnsF90()
437: !
438: if (fd_coloring) then
439: call MatFDColoringGetPerturbedColumnsF90(fdcoloring,indices,ierr)
440: print*,'Indices from GetPerturbedColumnsF90'
441: print*,indices
442: call MatFDColoringRestorePerturbedColumnsF90(fdcoloring,indices,ierr)
443: endif
444: return
445: end
447: ! ---------------------------------------------------------------------
448: !
449: ! ApplicationFunction - Computes nonlinear function, called by
450: ! the higher level routine FormFunction().
451: !
452: ! Input Parameter:
453: ! x - local vector data
454: !
455: ! Output Parameters:
456: ! f - local vector data, f(x)
457: ! ierr - error code
458: !
459: ! Notes:
460: ! This routine uses standard Fortran-style computations over a 2-dim array.
461: !
462: subroutine ApplicationFunction(x,f,ierr)
463: use petscsnes
464: implicit none
466: ! Common blocks:
467: PetscReal lambda
468: PetscInt mx,my
469: MatFDColoring fdcoloring
470: PetscBool fd_coloring
471: common /params/ lambda,mx,my,fdcoloring,fd_coloring
473: ! Input/output variables:
474: PetscScalar x(mx,my),f(mx,my)
475: PetscErrorCode ierr
477: ! Local variables:
478: PetscScalar two,one,hx,hy
479: PetscScalar hxdhy,hydhx,sc
480: PetscScalar u,uxx,uyy
481: PetscInt i,j
483: 0
484: one = 1.0
485: two = 2.0
486: hx = one/(mx-1)
487: hy = one/(my-1)
488: sc = hx*hy*lambda
489: hxdhy = hx/hy
490: hydhx = hy/hx
492: ! Compute function
494: do 20 j=1,my
495: do 10 i=1,mx
496: if (i .eq. 1 .or. j .eq. 1 .or. i .eq. mx .or. j .eq. my) then
497: f(i,j) = x(i,j)
498: else
499: u = x(i,j)
500: uxx = hydhx * (two*u - x(i-1,j) - x(i+1,j))
501: uyy = hxdhy * (two*u - x(i,j-1) - x(i,j+1))
502: f(i,j) = uxx + uyy - sc*exp(u)
503: endif
504: 10 continue
505: 20 continue
507: return
508: end
510: ! ---------------------------------------------------------------------
511: !
512: ! FormJacobian - Evaluates Jacobian matrix.
513: !
514: ! Input Parameters:
515: ! snes - the SNES context
516: ! x - input vector
517: ! dummy - optional user-defined context, as set by SNESSetJacobian()
518: ! (not used here)
519: !
520: ! Output Parameters:
521: ! jac - Jacobian matrix
522: ! jac_prec - optionally different preconditioning matrix (not used here)
523: ! flag - flag indicating matrix structure
524: !
525: ! Notes:
526: ! This routine serves as a wrapper for the lower-level routine
527: ! "ApplicationJacobian", where the actual computations are
528: ! done using the standard Fortran style of treating the local
529: ! vector data as a multidimensional array over the local mesh.
530: ! This routine merely accesses the local vector data via
531: ! VecGetArray() and VecRestoreArray().
532: !
533: subroutine FormJacobian(snes,X,jac,jac_prec,dummy,ierr)
534: use petscsnes
535: implicit none
537: ! Input/output variables:
538: SNES snes
539: Vec X
540: Mat jac,jac_prec
541: PetscErrorCode ierr
542: integer dummy
544: ! Common blocks:
545: PetscReal lambda
546: PetscInt mx,my
547: MatFDColoring fdcoloring
548: PetscBool fd_coloring
549: common /params/ lambda,mx,my,fdcoloring,fd_coloring
551: ! Declarations for use with local array:
552: PetscScalar lx_v(0:1)
553: PetscOffset lx_i
555: ! Get a pointer to vector data
557: call VecGetArrayRead(X,lx_v,lx_i,ierr)
559: ! Compute Jacobian entries
561: call ApplicationJacobian(lx_v(lx_i),jac,jac_prec,ierr)
563: ! Restore vector
565: call VecRestoreArrayRead(X,lx_v,lx_i,ierr)
567: ! Assemble matrix
569: call MatAssemblyBegin(jac_prec,MAT_FINAL_ASSEMBLY,ierr)
570: call MatAssemblyEnd(jac_prec,MAT_FINAL_ASSEMBLY,ierr)
573: return
574: end
576: ! ---------------------------------------------------------------------
577: !
578: ! ApplicationJacobian - Computes Jacobian matrix, called by
579: ! the higher level routine FormJacobian().
580: !
581: ! Input Parameters:
582: ! x - local vector data
583: !
584: ! Output Parameters:
585: ! jac - Jacobian matrix
586: ! jac_prec - optionally different preconditioning matrix (not used here)
587: ! ierr - error code
588: !
589: ! Notes:
590: ! This routine uses standard Fortran-style computations over a 2-dim array.
591: !
592: subroutine ApplicationJacobian(x,jac,jac_prec,ierr)
593: use petscsnes
594: implicit none
596: ! Common blocks:
597: PetscReal lambda
598: PetscInt mx,my
599: MatFDColoring fdcoloring
600: PetscBool fd_coloring
601: common /params/ lambda,mx,my,fdcoloring,fd_coloring
603: ! Input/output variables:
604: PetscScalar x(mx,my)
605: Mat jac,jac_prec
606: PetscErrorCode ierr
608: ! Local variables:
609: PetscInt i,j,row(1),col(5),i1,i5
610: PetscScalar two,one, hx,hy
611: PetscScalar hxdhy,hydhx,sc,v(5)
613: ! Set parameters
615: i1 = 1
616: i5 = 5
617: one = 1.0
618: two = 2.0
619: hx = one/(mx-1)
620: hy = one/(my-1)
621: sc = hx*hy
622: hxdhy = hx/hy
623: hydhx = hy/hx
625: ! Compute entries of the Jacobian matrix
626: ! - Here, we set all entries for a particular row at once.
627: ! - Note that MatSetValues() uses 0-based row and column numbers
628: ! in Fortran as well as in C.
630: do 20 j=1,my
631: row(1) = (j-1)*mx - 1
632: do 10 i=1,mx
633: row(1) = row(1) + 1
634: ! boundary points
635: if (i .eq. 1 .or. j .eq. 1 .or. i .eq. mx .or. j .eq. my) then
636: call MatSetValues(jac_prec,i1,row,i1,row,one,INSERT_VALUES,ierr)
637: ! interior grid points
638: else
639: v(1) = -hxdhy
640: v(2) = -hydhx
641: v(3) = two*(hydhx + hxdhy) - sc*lambda*exp(x(i,j))
642: v(4) = -hydhx
643: v(5) = -hxdhy
644: col(1) = row(1) - mx
645: col(2) = row(1) - 1
646: col(3) = row(1)
647: col(4) = row(1) + 1
648: col(5) = row(1) + mx
649: call MatSetValues(jac_prec,i1,row,i5,col,v,INSERT_VALUES,ierr)
650: endif
651: 10 continue
652: 20 continue
654: return
655: end