Actual source code: ex14f.F
petsc-3.6.4 2016-04-12
1: !
2: !
3: ! Solves a nonlinear system in parallel with a user-defined
4: ! Newton method that uses KSP to solve the linearized Newton sytems. This solver
5: ! is a very simplistic inexact Newton method. The intent of this code is to
6: ! demonstrate the repeated solution of linear sytems with the same nonzero pattern.
7: !
8: ! This is NOT the recommended approach for solving nonlinear problems with PETSc!
9: ! We urge users to employ the SNES component for solving nonlinear problems whenever
10: ! possible, as it offers many advantages over coding nonlinear solvers independently.
11: !
12: ! We solve the Bratu (SFI - solid fuel ignition) problem in a 2D rectangular
13: ! domain, using distributed arrays (DMDAs) to partition the parallel grid.
14: !
15: ! The command line options include:
16: ! -par <parameter>, where <parameter> indicates the problem's nonlinearity
17: ! problem SFI: <parameter> = Bratu parameter (0 <= par <= 6.81)
18: ! -mx <xg>, where <xg> = number of grid points in the x-direction
19: ! -my <yg>, where <yg> = number of grid points in the y-direction
20: ! -Nx <npx>, where <npx> = number of processors in the x-direction
21: ! -Ny <npy>, where <npy> = number of processors in the y-direction
22: ! -mf use matrix free for matrix vector product
23: !
24: !/*T
25: ! Concepts: KSP^writing a user-defined nonlinear solver
26: ! Concepts: DMDA^using distributed arrays
27: ! Processors: n
28: !T*/
29: ! ------------------------------------------------------------------------
30: !
31: ! Solid Fuel Ignition (SFI) problem. This problem is modeled by
32: ! the partial differential equation
33: !
34: ! -Laplacian u - lambda*exp(u) = 0, 0 < x,y < 1,
35: !
36: ! with boundary conditions
37: !
38: ! u = 0 for x = 0, x = 1, y = 0, y = 1.
39: !
40: ! A finite difference approximation with the usual 5-point stencil
41: ! is used to discretize the boundary value problem to obtain a nonlinear
42: ! system of equations.
43: !
44: ! The SNES version of this problem is: snes/examples/tutorials/ex5f.F
45: !
46: ! -------------------------------------------------------------------------
48: program main
49: implicit none
51: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
52: ! Include files
53: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
54: !
55: ! petscsys.h - base PETSc routines petscvec.h - vectors
56: ! petscmat.h - matrices
57: ! petscis.h - index sets petscksp.h - Krylov subspace methods
58: ! petscviewer.h - viewers petscpc.h - preconditioners
60: #include <petsc/finclude/petscsys.h>
61: #include <petsc/finclude/petscis.h>
62: #include <petsc/finclude/petscvec.h>
63: #include <petsc/finclude/petscmat.h>
64: #include <petsc/finclude/petscpc.h>
65: #include <petsc/finclude/petscksp.h>
66: #include <petsc/finclude/petscdm.h>
67: #include <petsc/finclude/petscdmda.h>
69: MPI_Comm comm
70: Vec X,Y,F,localX
71: Mat J,B
72: DM da
73: KSP ksp
75: PetscInt Nx,Ny,N,mx,my,ifive,ithree
76: PetscBool flg,nooutput,usemf
77: common /mycommon/ mx,my,B,localX,da
78: !
79: !
80: ! This is the routine to use for matrix-free approach
81: !
82: external mymult
84: ! --------------- Data to define nonlinear solver --------------
85: PetscReal rtol,ttol
86: PetscReal fnorm,ynorm,xnorm
87: PetscInt max_nonlin_its,one
88: PetscInt lin_its
89: PetscInt i,m
90: PetscScalar mone
91: PetscErrorCode ierr
93: mone = -1.d0
94: rtol = 1.d-8
95: max_nonlin_its = 10
96: one = 1
97: ifive = 5
98: ithree = 3
100: call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
101: comm = PETSC_COMM_WORLD
103: ! Initialize problem parameters
105: !
106: mx = 4
107: my = 4
108: call PetscOptionsGetInt(PETSC_NULL_CHARACTER,'-mx',mx,flg,ierr)
109: call PetscOptionsGetInt(PETSC_NULL_CHARACTER,'-my',my,flg,ierr)
110: N = mx*my
112: nooutput = .false.
113: call PetscOptionsHasName(PETSC_NULL_CHARACTER,'-no_output', &
114: & nooutput,ierr)
116: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
117: ! Create linear solver context
118: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
120: call KSPCreate(comm,ksp,ierr)
122: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
123: ! Create vector data structures
124: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
126: !
127: ! Create distributed array (DMDA) to manage parallel grid and vectors
128: !
129: Nx = PETSC_DECIDE
130: Ny = PETSC_DECIDE
131: call PetscOptionsGetInt(PETSC_NULL_CHARACTER,'-Nx',Nx,flg,ierr)
132: call PetscOptionsGetInt(PETSC_NULL_CHARACTER,'-Ny',Ny,flg,ierr)
133: call DMDACreate2d(comm,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE, &
134: & DMDA_STENCIL_STAR,mx,my,Nx,Ny,one,one, &
135: & PETSC_NULL_INTEGER,PETSC_NULL_INTEGER,da,ierr)
137: !
138: ! Extract global and local vectors from DMDA then duplicate for remaining
139: ! vectors that are the same types
140: !
141: call DMCreateGlobalVector(da,X,ierr)
142: call DMCreateLocalVector(da,localX,ierr)
143: call VecDuplicate(X,F,ierr)
144: call VecDuplicate(X,Y,ierr)
147: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
148: ! Create matrix data structure for Jacobian
149: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
150: !
151: ! Note: For the parallel case, vectors and matrices MUST be partitioned
152: ! accordingly. When using distributed arrays (DMDAs) to create vectors,
153: ! the DMDAs determine the problem partitioning. We must explicitly
154: ! specify the local matrix dimensions upon its creation for compatibility
155: ! with the vector distribution.
156: !
157: ! Note: Here we only approximately preallocate storage space for the
158: ! Jacobian. See the users manual for a discussion of better techniques
159: ! for preallocating matrix memory.
160: !
161: call VecGetLocalSize(X,m,ierr)
162: call MatCreateAIJ(comm,m,m,N,N,ifive,PETSC_NULL_INTEGER,ithree, &
163: & PETSC_NULL_INTEGER,B,ierr)
165: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
166: ! if usemf is on then matrix vector product is done via matrix free
167: ! approach. Note this is just an example, and not realistic because
168: ! we still use the actual formed matrix, but in reality one would
169: ! provide their own subroutine that would directly do the matrix
170: ! vector product and not call MatMult()
171: ! Note: we put B into a common block so it will be visible to the
172: ! mymult() routine
173: usemf = .false.
174: call PetscOptionsHasName(PETSC_NULL_CHARACTER,'-mf',usemf,ierr)
175: if (usemf) then
176: call MatCreateShell(comm,m,m,N,N,PETSC_NULL_INTEGER,J,ierr)
177: call MatShellSetOperation(J,MATOP_MULT,mymult,ierr)
178: else
179: ! If not doing matrix free then matrix operator, J, and matrix used
180: ! to construct preconditioner, B, are the same
181: J = B
182: endif
184: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
185: ! Customize linear solver set runtime options
186: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
187: !
188: ! Set runtime options (e.g., -ksp_monitor -ksp_rtol <rtol> -ksp_type <type>)
189: !
190: call KSPSetFromOptions(ksp,ierr)
192: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
193: ! Evaluate initial guess
194: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
196: call FormInitialGuess(X,ierr)
197: call ComputeFunction(X,F,ierr)
198: call VecNorm(F,NORM_2,fnorm,ierr)
199: ttol = fnorm*rtol
200: if (.not. nooutput) then
201: print*, 'Initial function norm ',fnorm
202: endif
204: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
205: ! Solve nonlinear system with a user-defined method
206: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
208: ! This solver is a very simplistic inexact Newton method, with no
209: ! no damping strategies or bells and whistles. The intent of this code
210: ! is merely to demonstrate the repeated solution with KSP of linear
211: ! sytems with the same nonzero structure.
212: !
213: ! This is NOT the recommended approach for solving nonlinear problems
214: ! with PETSc! We urge users to employ the SNES component for solving
215: ! nonlinear problems whenever possible with application codes, as it
216: ! offers many advantages over coding nonlinear solvers independently.
218: do 10 i=0,max_nonlin_its
220: ! Compute the Jacobian matrix. See the comments in this routine for
221: ! important information about setting the flag mat_flag.
223: call ComputeJacobian(X,B,ierr)
225: ! Solve J Y = F, where J is the Jacobian matrix.
226: ! - First, set the KSP linear operators. Here the matrix that
227: ! defines the linear system also serves as the preconditioning
228: ! matrix.
229: ! - Then solve the Newton system.
231: call KSPSetOperators(ksp,J,B,ierr)
232: call KSPSolve(ksp,F,Y,ierr)
234: ! Compute updated iterate
236: call VecNorm(Y,NORM_2,ynorm,ierr)
237: call VecAYPX(Y,mone,X,ierr)
238: call VecCopy(Y,X,ierr)
239: call VecNorm(X,NORM_2,xnorm,ierr)
240: call KSPGetIterationNumber(ksp,lin_its,ierr)
241: if (.not. nooutput) then
242: print*,'linear solve iterations = ',lin_its,' xnorm = ', &
243: & xnorm,' ynorm = ',ynorm
244: endif
246: ! Evaluate nonlinear function at new location
248: call ComputeFunction(X,F,ierr)
249: call VecNorm(F,NORM_2,fnorm,ierr)
250: if (.not. nooutput) then
251: print*, 'Iteration ',i+1,' function norm',fnorm
252: endif
254: ! Test for convergence
256: if (fnorm .le. ttol) then
257: if (.not. nooutput) then
258: print*,'Converged: function norm ',fnorm,' tolerance ',ttol
259: endif
260: goto 20
261: endif
262: 10 continue
263: 20 continue
265: write(6,100) i+1
266: 100 format('Number of SNES iterations =',I2)
268: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
269: ! Free work space. All PETSc objects should be destroyed when they
270: ! are no longer needed.
271: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
273: call MatDestroy(B,ierr)
274: if (usemf) then
275: call MatDestroy(J,ierr)
276: endif
277: call VecDestroy(localX,ierr)
278: call VecDestroy(X,ierr)
279: call VecDestroy(Y,ierr)
280: call VecDestroy(F,ierr)
281: call KSPDestroy(ksp,ierr)
282: call DMDestroy(da,ierr)
283: call PetscFinalize(ierr)
284: end
286: ! -------------------------------------------------------------------
287: !
288: ! FormInitialGuess - Forms initial approximation.
289: !
290: ! Input Parameters:
291: ! X - vector
292: !
293: ! Output Parameter:
294: ! X - vector
295: !
296: subroutine FormInitialGuess(X,ierr)
297: implicit none
299: ! petscsys.h - base PETSc routines petscvec.h - vectors
300: ! petscmat.h - matrices
301: ! petscis.h - index sets petscksp.h - Krylov subspace methods
302: ! petscviewer.h - viewers petscpc.h - preconditioners
304: #include <petsc/finclude/petscsys.h>
305: #include <petsc/finclude/petscis.h>
306: #include <petsc/finclude/petscvec.h>
307: #include <petsc/finclude/petscmat.h>
308: #include <petsc/finclude/petscpc.h>
309: #include <petsc/finclude/petscksp.h>
310: #include <petsc/finclude/petscdm.h>
311: #include <petsc/finclude/petscdmda.h>
312: PetscErrorCode ierr
313: PetscOffset idx
314: Vec X,localX
315: PetscInt i,j,row,mx
316: PetscInt my, xs,ys,xm
317: PetscInt ym
318: PetscReal one,lambda,temp1,temp,hx,hy
319: PetscScalar xx(1)
320: DM da
321: Mat B
322: common /mycommon/ mx,my,B,localX,da
324: one = 1.d0
325: lambda = 6.d0
326: hx = one/(mx-1)
327: hy = one/(my-1)
328: temp1 = lambda/(lambda + one)
330: ! Get a pointer to vector data.
331: ! - VecGetArray() returns a pointer to the data array.
332: ! - You MUST call VecRestoreArray() when you no longer need access to
333: ! the array.
334: call VecGetArray(X,xx,idx,ierr)
336: ! Get local grid boundaries (for 2-dimensional DMDA):
337: ! xs, ys - starting grid indices (no ghost points)
338: ! xm, ym - widths of local grid (no ghost points)
340: call DMDAGetCorners(da,xs,ys,PETSC_NULL_INTEGER,xm,ym, &
341: & PETSC_NULL_INTEGER,ierr)
343: ! Compute initial guess over the locally owned part of the grid
345: do 30 j=ys,ys+ym-1
346: temp = (min(j,my-j-1))*hy
347: do 40 i=xs,xs+xm-1
348: row = i - xs + (j - ys)*xm + 1
349: if (i .eq. 0 .or. j .eq. 0 .or. i .eq. mx-1 .or. &
350: & j .eq. my-1) then
351: xx(idx+row) = 0.d0
352: continue
353: endif
354: xx(idx+row) = temp1*sqrt(min((min(i,mx-i-1))*hx,temp))
355: 40 continue
356: 30 continue
358: ! Restore vector
360: call VecRestoreArray(X,xx,idx,ierr)
361: return
362: end
364: ! -------------------------------------------------------------------
365: !
366: ! ComputeFunction - Evaluates nonlinear function, F(x).
367: !
368: ! Input Parameters:
369: !. X - input vector
370: !
371: ! Output Parameter:
372: !. F - function vector
373: !
374: subroutine ComputeFunction(X,F,ierr)
375: implicit none
377: ! petscsys.h - base PETSc routines petscvec.h - vectors
378: ! petscmat.h - matrices
379: ! petscis.h - index sets petscksp.h - Krylov subspace methods
380: ! petscviewer.h - viewers petscpc.h - preconditioners
382: #include <petsc/finclude/petscsys.h>
383: #include <petsc/finclude/petscis.h>
384: #include <petsc/finclude/petscvec.h>
385: #include <petsc/finclude/petscmat.h>
386: #include <petsc/finclude/petscpc.h>
387: #include <petsc/finclude/petscksp.h>
388: #include <petsc/finclude/petscdm.h>
389: #include <petsc/finclude/petscdmda.h>
391: Vec X,F,localX
392: PetscInt gys,gxm,gym
393: PetscOffset idx,idf
394: PetscErrorCode ierr
395: PetscInt i,j,row,mx,my,xs,ys,xm,ym,gxs
396: PetscInt rowf
397: PetscReal two,one,lambda,hx
398: PetscReal hy,hxdhy,hydhx,sc
399: PetscScalar u,uxx,uyy,xx(1),ff(1)
400: DM da
401: Mat B
402: common /mycommon/ mx,my,B,localX,da
404: two = 2.d0
405: one = 1.d0
406: lambda = 6.d0
408: hx = one/(mx-1)
409: hy = one/(my-1)
410: sc = hx*hy*lambda
411: hxdhy = hx/hy
412: hydhx = hy/hx
414: ! Scatter ghost points to local vector, using the 2-step process
415: ! DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
416: ! By placing code between these two statements, computations can be
417: ! done while messages are in transition.
418: !
419: call DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX,ierr)
420: call DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX,ierr)
422: ! Get pointers to vector data
424: call VecGetArray(localX,xx,idx,ierr)
425: call VecGetArray(F,ff,idf,ierr)
427: ! Get local grid boundaries
429: call DMDAGetCorners(da,xs,ys,PETSC_NULL_INTEGER,xm,ym, &
430: & PETSC_NULL_INTEGER,ierr)
431: call DMDAGetGhostCorners(da,gxs,gys,PETSC_NULL_INTEGER,gxm,gym, &
432: & PETSC_NULL_INTEGER,ierr)
434: ! Compute function over the locally owned part of the grid
435: rowf = 0
436: do 50 j=ys,ys+ym-1
438: row = (j - gys)*gxm + xs - gxs
439: do 60 i=xs,xs+xm-1
440: row = row + 1
441: rowf = rowf + 1
443: if (i .eq. 0 .or. j .eq. 0 .or. i .eq. mx-1 .or. &
444: & j .eq. my-1) then
445: ff(idf+rowf) = xx(idx+row)
446: goto 60
447: endif
448: u = xx(idx+row)
449: uxx = (two*u - xx(idx+row-1) - xx(idx+row+1))*hydhx
450: uyy = (two*u - xx(idx+row-gxm) - xx(idx+row+gxm))*hxdhy
451: ff(idf+rowf) = uxx + uyy - sc*exp(u)
452: 60 continue
453: 50 continue
455: ! Restore vectors
457: call VecRestoreArray(localX,xx,idx,ierr)
458: call VecRestoreArray(F,ff,idf,ierr)
459: return
460: end
462: ! -------------------------------------------------------------------
463: !
464: ! ComputeJacobian - Evaluates Jacobian matrix.
465: !
466: ! Input Parameters:
467: ! x - input vector
468: !
469: ! Output Parameters:
470: ! jac - Jacobian matrix
471: ! flag - flag indicating matrix structure
472: !
473: ! Notes:
474: ! Due to grid point reordering with DMDAs, we must always work
475: ! with the local grid points, and then transform them to the new
476: ! global numbering with the 'ltog' mapping
477: ! We cannot work directly with the global numbers for the original
478: ! uniprocessor grid!
479: !
480: subroutine ComputeJacobian(X,jac,ierr)
481: implicit none
483: ! petscsys.h - base PETSc routines petscvec.h - vectors
484: ! petscmat.h - matrices
485: ! petscis.h - index sets petscksp.h - Krylov subspace methods
486: ! petscviewer.h - viewers petscpc.h - preconditioners
488: #include <petsc/finclude/petscsys.h>
489: #include <petsc/finclude/petscis.h>
490: #include <petsc/finclude/petscvec.h>
491: #include <petsc/finclude/petscmat.h>
492: #include <petsc/finclude/petscpc.h>
493: #include <petsc/finclude/petscksp.h>
494: #include <petsc/finclude/petscdm.h>
495: #include <petsc/finclude/petscdmda.h>
497: Vec X
498: Mat jac
499: Vec localX
500: DM da
501: PetscInt ltog(1)
502: PetscOffset idltog,idx
503: PetscErrorCode ierr
504: PetscInt xs,ys,xm,ym
505: PetscInt gxs,gys,gxm,gym
506: PetscInt grow(1),i,j
507: PetscInt row,mx,my,ione
508: PetscInt col(5),ifive
509: PetscScalar two,one,lambda
510: PetscScalar v(5),hx,hy,hxdhy
511: PetscScalar hydhx,sc,xx(1)
512: Mat B
513: ISLocalToGlobalMapping ltogm
514: common /mycommon/ mx,my,B,localX,da
516: ione = 1
517: ifive = 5
518: one = 1.d0
519: two = 2.d0
520: hx = one/(mx-1)
521: hy = one/(my-1)
522: sc = hx*hy
523: hxdhy = hx/hy
524: hydhx = hy/hx
525: lambda = 6.d0
527: ! Scatter ghost points to local vector, using the 2-step process
528: ! DMGlobalToLocalBegin(), DMGlobalToLocalEnd().
529: ! By placing code between these two statements, computations can be
530: ! done while messages are in transition.
532: call DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX,ierr)
533: call DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX,ierr)
535: ! Get pointer to vector data
537: call VecGetArray(localX,xx,idx,ierr)
539: ! Get local grid boundaries
541: call DMDAGetCorners(da,xs,ys,PETSC_NULL_INTEGER,xm,ym, &
542: & PETSC_NULL_INTEGER,ierr)
543: call DMDAGetGhostCorners(da,gxs,gys,PETSC_NULL_INTEGER,gxm,gym, &
544: & PETSC_NULL_INTEGER,ierr)
546: ! Get the global node numbers for all local nodes, including ghost points
548: call DMGetLocalToGlobalMapping(da,ltogm,ierr)
549: call ISLocalToGlobalMappingGetIndices(ltogm,ltog,idltog,ierr)
551: ! Compute entries for the locally owned part of the Jacobian.
552: ! - Currently, all PETSc parallel matrix formats are partitioned by
553: ! contiguous chunks of rows across the processors. The 'grow'
554: ! parameter computed below specifies the global row number
555: ! corresponding to each local grid point.
556: ! - Each processor needs to insert only elements that it owns
557: ! locally (but any non-local elements will be sent to the
558: ! appropriate processor during matrix assembly).
559: ! - Always specify global row and columns of matrix entries.
560: ! - Here, we set all entries for a particular row at once.
562: do 10 j=ys,ys+ym-1
563: row = (j - gys)*gxm + xs - gxs
564: do 20 i=xs,xs+xm-1
565: row = row + 1
566: grow(1) = ltog(idltog+row)
567: if (i .eq. 0 .or. j .eq. 0 .or. i .eq. (mx-1) .or. &
568: & j .eq. (my-1)) then
569: call MatSetValues(jac,ione,grow,ione,grow,one, &
570: & INSERT_VALUES,ierr)
571: go to 20
572: endif
573: v(1) = -hxdhy
574: col(1) = ltog(idltog+row - gxm)
575: v(2) = -hydhx
576: col(2) = ltog(idltog+row - 1)
577: v(3) = two*(hydhx + hxdhy) - sc*lambda*exp(xx(idx+row))
578: col(3) = grow(1)
579: v(4) = -hydhx
580: col(4) = ltog(idltog+row + 1)
581: v(5) = -hxdhy
582: col(5) = ltog(idltog+row + gxm)
583: call MatSetValues(jac,ione,grow,ifive,col,v,INSERT_VALUES, &
584: & ierr)
585: 20 continue
586: 10 continue
588: call ISLocalToGlobalMappingRestoreIndices(ltogm,ltog,idltog,ierr)
590: ! Assemble matrix, using the 2-step process:
591: ! MatAssemblyBegin(), MatAssemblyEnd().
592: ! By placing code between these two statements, computations can be
593: ! done while messages are in transition.
595: call MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY,ierr)
596: call VecRestoreArray(localX,xx,idx,ierr)
597: call MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY,ierr)
598: return
599: end
602: ! -------------------------------------------------------------------
603: !
604: ! MyMult - user provided matrix multiply
605: !
606: ! Input Parameters:
607: !. X - input vector
608: !
609: ! Output Parameter:
610: !. F - function vector
611: !
612: subroutine MyMult(J,X,F,ierr)
613: implicit none
614: Mat J,B
615: Vec X,F
616: PetscErrorCode ierr
617: PetscInt mx,my
618: DM da
619: Vec localX
621: common /mycommon/ mx,my,B,localX,da
622: !
623: ! Here we use the actual formed matrix B; users would
624: ! instead write their own matrix vector product routine
625: !
626: call MatMult(B,X,F,ierr)
627: return
628: end