Actual source code: ex1f.F90
petsc-3.9.0 2018-04-07
1: !
2: !
3: ! Description: Uses the Newton method to solve a two-variable system.
4: !
5: !!/*T
6: ! Concepts: SNES^basic uniprocessor example
7: ! Processors: 1
8: !T*/
11: !
12: ! -----------------------------------------------------------------------
14: program main
15: #include <petsc/finclude/petsc.h>
16: use petsc
17: implicit none
20: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
21: ! Variable declarations
22: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
23: !
24: ! Variables:
25: ! snes - nonlinear solver
26: ! ksp - linear solver
27: ! pc - preconditioner context
28: ! ksp - Krylov subspace method context
29: ! x, r - solution, residual vectors
30: ! J - Jacobian matrix
31: ! its - iterations for convergence
32: !
33: SNES snes
34: PC pc
35: KSP ksp
36: Vec x,r
37: Mat J
38: SNESLineSearch linesearch
39: PetscErrorCode ierr
40: PetscInt its,i2,i20
41: PetscMPIInt size,rank
42: PetscScalar pfive
43: PetscReal tol
44: PetscBool setls
46: ! Note: Any user-defined Fortran routines (such as FormJacobian)
47: ! MUST be declared as external.
49: external FormFunction, FormJacobian, MyLineSearch
51: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
52: ! Macro definitions
53: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
54: !
55: ! Macros to make clearer the process of setting values in vectors and
56: ! getting values from vectors. These vectors are used in the routines
57: ! FormFunction() and FormJacobian().
58: ! - The element lx_a(ib) is element ib in the vector x
59: !
60: #define lx_a(ib) lx_v(lx_i + (ib))
61: #define lf_a(ib) lf_v(lf_i + (ib))
62: !
63: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
64: ! Beginning of program
65: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
67: call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
68: if (ierr .ne. 0) then
69: print*,'Unable to initialize PETSc'
70: stop
71: endif
72: call MPI_Comm_size(PETSC_COMM_WORLD,size,ierr)
73: call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
74: if (size .ne. 1) then; SETERRA(PETSC_COMM_SELF,1,'Uniprocessor example'); endif
76: i2 = 2
77: i20 = 20
78: ! - - - - - - - - - -- - - - - - - - - - - - - - - - - - - - - - - - - -
79: ! Create nonlinear solver context
80: ! - - - - - - - - - -- - - - - - - - - - - - - - - - - - - - - - - - - -
82: call SNESCreate(PETSC_COMM_WORLD,snes,ierr)
84: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
85: ! Create matrix and vector data structures; set corresponding routines
86: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
88: ! Create vectors for solution and nonlinear function
90: call VecCreateSeq(PETSC_COMM_SELF,i2,x,ierr)
91: call VecDuplicate(x,r,ierr)
93: ! Create Jacobian matrix data structure
95: call MatCreate(PETSC_COMM_SELF,J,ierr)
96: call MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,i2,i2,ierr)
97: call MatSetFromOptions(J,ierr)
98: call MatSetUp(J,ierr)
100: ! Set function evaluation routine and vector
102: call SNESSetFunction(snes,r,FormFunction,0,ierr)
104: ! Set Jacobian matrix data structure and Jacobian evaluation routine
106: call SNESSetJacobian(snes,J,J,FormJacobian,0,ierr)
108: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
109: ! Customize nonlinear solver; set runtime options
110: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
112: ! Set linear solver defaults for this problem. By extracting the
113: ! KSP, KSP, and PC contexts from the SNES context, we can then
114: ! directly call any KSP, KSP, and PC routines to set various options.
116: call SNESGetKSP(snes,ksp,ierr)
117: call KSPGetPC(ksp,pc,ierr)
118: call PCSetType(pc,PCNONE,ierr)
119: tol = 1.e-4
120: call KSPSetTolerances(ksp,tol,PETSC_DEFAULT_REAL, &
121: & PETSC_DEFAULT_REAL,i20,ierr)
123: ! Set SNES/KSP/KSP/PC runtime options, e.g.,
124: ! -snes_view -snes_monitor -ksp_type <ksp> -pc_type <pc>
125: ! These options will override those specified above as long as
126: ! SNESSetFromOptions() is called _after_ any other customization
127: ! routines.
130: call SNESSetFromOptions(snes,ierr)
132: call PetscOptionsHasName(PETSC_NULL_OPTIONS,PETSC_NULL_CHARACTER, &
133: & '-setls',setls,ierr)
135: if (setls) then
136: call SNESGetLineSearch(snes, linesearch, ierr)
137: call SNESLineSearchSetType(linesearch, 'shell', ierr)
138: call SNESLineSearchShellSetUserFunc(linesearch, MyLineSearch, &
139: & 0, ierr)
140: endif
142: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
143: ! Evaluate initial guess; then solve nonlinear system
144: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
146: ! Note: The user should initialize the vector, x, with the initial guess
147: ! for the nonlinear solver prior to calling SNESSolve(). In particular,
148: ! to employ an initial guess of zero, the user should explicitly set
149: ! this vector to zero by calling VecSet().
151: pfive = 0.5
152: call VecSet(x,pfive,ierr)
153: call SNESSolve(snes,PETSC_NULL_VEC,x,ierr)
154: call SNESGetIterationNumber(snes,its,ierr);
155: if (rank .eq. 0) then
156: write(6,100) its
157: endif
158: 100 format('Number of SNES iterations = ',i5)
160: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
161: ! Free work space. All PETSc objects should be destroyed when they
162: ! are no longer needed.
163: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
165: call VecDestroy(x,ierr)
166: call VecDestroy(r,ierr)
167: call MatDestroy(J,ierr)
168: call SNESDestroy(snes,ierr)
169: call PetscFinalize(ierr)
170: end
171: !
172: ! ------------------------------------------------------------------------
173: !
174: ! FormFunction - Evaluates nonlinear function, F(x).
175: !
176: ! Input Parameters:
177: ! snes - the SNES context
178: ! x - input vector
179: ! dummy - optional user-defined context (not used here)
180: !
181: ! Output Parameter:
182: ! f - function vector
183: !
184: subroutine FormFunction(snes,x,f,dummy,ierr)
185: use petscsnes
186: implicit none
188: SNES snes
189: Vec x,f
190: PetscErrorCode ierr
191: integer dummy(*)
193: ! Declarations for use with local arrays
195: PetscScalar lx_v(2),lf_v(2)
196: PetscOffset lx_i,lf_i
198: ! Get pointers to vector data.
199: ! - For default PETSc vectors, VecGetArray() returns a pointer to
200: ! the data array. Otherwise, the routine is implementation dependent.
201: ! - You MUST call VecRestoreArray() when you no longer need access to
202: ! the array.
203: ! - Note that the Fortran interface to VecGetArray() differs from the
204: ! C version. See the Fortran chapter of the users manual for details.
206: call VecGetArrayRead(x,lx_v,lx_i,ierr)
207: call VecGetArray(f,lf_v,lf_i,ierr)
209: ! Compute function
211: lf_a(1) = lx_a(1)*lx_a(1) &
212: & + lx_a(1)*lx_a(2) - 3.0
213: lf_a(2) = lx_a(1)*lx_a(2) &
214: & + lx_a(2)*lx_a(2) - 6.0
216: ! Restore vectors
218: call VecRestoreArrayRead(x,lx_v,lx_i,ierr)
219: call VecRestoreArray(f,lf_v,lf_i,ierr)
221: return
222: end
224: ! ---------------------------------------------------------------------
225: !
226: ! FormJacobian - Evaluates Jacobian matrix.
227: !
228: ! Input Parameters:
229: ! snes - the SNES context
230: ! x - input vector
231: ! dummy - optional user-defined context (not used here)
232: !
233: ! Output Parameters:
234: ! A - Jacobian matrix
235: ! B - optionally different preconditioning matrix
236: ! flag - flag indicating matrix structure
237: !
238: subroutine FormJacobian(snes,X,jac,B,dummy,ierr)
239: use petscsnes
240: implicit none
242: SNES snes
243: Vec X
244: Mat jac,B
245: PetscScalar A(4)
246: PetscErrorCode ierr
247: PetscInt idx(2),i2
248: integer dummy(*)
250: ! Declarations for use with local arrays
252: PetscScalar lx_v(2)
253: PetscOffset lx_i
255: ! Get pointer to vector data
257: i2 = 2
258: call VecGetArrayRead(x,lx_v,lx_i,ierr)
260: ! Compute Jacobian entries and insert into matrix.
261: ! - Since this is such a small problem, we set all entries for
262: ! the matrix at once.
263: ! - Note that MatSetValues() uses 0-based row and column numbers
264: ! in Fortran as well as in C (as set here in the array idx).
266: idx(1) = 0
267: idx(2) = 1
268: A(1) = 2.0*lx_a(1) + lx_a(2)
269: A(2) = lx_a(1)
270: A(3) = lx_a(2)
271: A(4) = lx_a(1) + 2.0*lx_a(2)
272: call MatSetValues(B,i2,idx,i2,idx,A,INSERT_VALUES,ierr)
274: ! Restore vector
276: call VecRestoreArrayRead(x,lx_v,lx_i,ierr)
278: ! Assemble matrix
280: call MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY,ierr)
281: call MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY,ierr)
282: if (B .ne. jac) then
283: call MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY,ierr)
284: call MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY,ierr)
285: endif
287: return
288: end
291: subroutine MyLineSearch(linesearch, lctx, ierr)
292: use petscsnes
293: implicit none
295: SNESLineSearch linesearch
296: SNES snes
297: integer lctx
298: Vec x, f,g, y, w
299: PetscReal ynorm,gnorm,xnorm
300: PetscBool flag
301: PetscErrorCode ierr
303: PetscScalar mone
305: mone = -1.0
306: call SNESLineSearchGetSNES(linesearch, snes, ierr)
307: call SNESLineSearchGetVecs(linesearch, x, f, y, w, g, ierr)
308: call VecNorm(y,NORM_2,ynorm,ierr)
309: call VecAXPY(x,mone,y,ierr)
310: call SNESComputeFunction(snes,x,f,ierr)
311: call VecNorm(f,NORM_2,gnorm,ierr)
312: call VecNorm(x,NORM_2,xnorm,ierr)
313: call VecNorm(y,NORM_2,ynorm,ierr)
314: call SNESLineSearchSetNorms(linesearch, xnorm, gnorm, ynorm, &
315: & ierr)
316: flag = PETSC_FALSE
317: return
318: end
320: !/*TEST
321: !
322: ! test:
323: ! args: -ksp_gmres_cgs_refinement_type refine_always -snes_monitor_short
324: ! requires: !single
325: !
326: !TEST*/