1: !
2: !
3: !/*T
4: ! Concepts: KSP^basic sequential example
5: ! Concepts: KSP^Laplacian, 2d
6: ! Concepts: Laplacian, 2d
7: ! Processors: 1
8: !T*/
9: ! -----------------------------------------------------------------------
11: program main
12: implicit none
14: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
15: ! Include files
16: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
17: !
18: ! The following include statements are required for KSP Fortran programs:
19: ! petscsys.h - base PETSc routines
20: ! petscvec.h - vectors
21: ! petscmat.h - matrices
22: ! petscksp.h - Krylov subspace methods
23: ! petscpc.h - preconditioners
24: !
25: #include <finclude/petscsys.h>
26: #include <finclude/petscvec.h>
27: #include <finclude/petscmat.h>
28: #include <finclude/petscksp.h>
29: #include <finclude/petscpc.h>
31: ! User-defined context that contains all the data structures used
32: ! in the linear solution process.
34: ! Vec x,b /* solution vector, right hand side vector and work vector */
35: ! Mat A /* sparse matrix */
36: ! KSP ksp /* linear solver context */
37: ! int m,n /* grid dimensions */
38: !
39: ! Since we cannot store Scalars and integers in the same context,
40: ! we store the integers/pointers in the user-defined context, and
41: ! the scalar values are carried in the common block.
42: ! The scalar values in this simplistic example could easily
43: ! be recalculated in each routine, where they are needed.
44: !
45: ! Scalar hx2,hy2 /* 1/(m+1)*(m+1) and 1/(n+1)*(n+1) */
47: ! Note: Any user-defined Fortran routines MUST be declared as external.
49: external UserInitializeLinearSolver
50: external UserFinalizeLinearSolver
51: external UserDoLinearSolver
53: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
54: ! Variable declarations
55: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
57: PetscScalar hx,hy,x,y
58: PetscFortranAddr userctx(6)
59: PetscErrorCode ierr
60: PetscInt m,n,t,tmax,i,j
61: PetscBool flg
62: PetscMPIInt size,rank
63: PetscReal enorm
64: PetscScalar cnorm
65: PetscScalar,ALLOCATABLE :: userx(:,:)
66: PetscScalar,ALLOCATABLE :: userb(:,:)
67: PetscScalar,ALLOCATABLE :: solution(:,:)
68: PetscScalar,ALLOCATABLE :: rho(:,:)
70: PetscReal hx2,hy2
71: common /param/ hx2,hy2
73: tmax = 2
74: m = 6
75: n = 7
77: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
78: ! Beginning of program
79: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
81: call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
82: call MPI_Comm_size(PETSC_COMM_WORLD,size,ierr)
83: if (size .ne. 1) then
84: call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
85: if (rank .eq. 0) then
86: write(6,*) 'This is a uniprocessor example only!'
87: endif
88: SETERRQ(PETSC_COMM_WORLD,1,' ',ierr)
89: endif
91: ! The next two lines are for testing only; these allow the user to
92: ! decide the grid size at runtime.
94: call PetscOptionsGetInt(PETSC_NULL_CHARACTER,'-m',m,flg,ierr)
95: call PetscOptionsGetInt(PETSC_NULL_CHARACTER,'-n',n,flg,ierr)
97: ! Create the empty sparse matrix and linear solver data structures
99: call UserInitializeLinearSolver(m,n,userctx,ierr)
101: ! Allocate arrays to hold the solution to the linear system. This
102: ! approach is not normally done in PETSc programs, but in this case,
103: ! since we are calling these routines from a non-PETSc program, we
104: ! would like to reuse the data structures from another code. So in
105: ! the context of a larger application these would be provided by
106: ! other (non-PETSc) parts of the application code.
108: ALLOCATE (userx(m,n),userb(m,n),solution(m,n))
110: ! Allocate an array to hold the coefficients in the elliptic operator
112: ALLOCATE (rho(m,n))
114: ! Fill up the array rho[] with the function rho(x,y) = x; fill the
115: ! right-hand-side b[] and the solution with a known problem for testing.
117: hx = 1.0/(m+1)
118: hy = 1.0/(n+1)
119: y = hy
120: do 20 j=1,n
121: x = hx
122: do 10 i=1,m
123: rho(i,j) = x
124: solution(i,j) = sin(2.*PETSC_PI*x)*sin(2.*PETSC_PI*y)
125: userb(i,j) = -2.*PETSC_PI*cos(2.*PETSC_PI*x)* &
126: & sin(2.*PETSC_PI*y) + &
127: & 8*PETSC_PI*PETSC_PI*x* &
128: & sin(2.*PETSC_PI*x)*sin(2.*PETSC_PI*y)
129: x = x + hx
130: 10 continue
131: y = y + hy
132: 20 continue
134: ! Loop over a bunch of timesteps, setting up and solver the linear
135: ! system for each time-step.
136: ! Note that this loop is somewhat artificial. It is intended to
137: ! demonstrate how one may reuse the linear solvers in each time-step.
139: do 100 t=1,tmax
140: call UserDoLinearSolver(rho,userctx,userb,userx,ierr)
142: ! Compute error: Note that this could (and usually should) all be done
143: ! using the PETSc vector operations. Here we demonstrate using more
144: ! standard programming practices to show how they may be mixed with
145: ! PETSc.
146: cnorm = 0.0
147: do 90 j=1,n
148: do 80 i=1,m
149: cnorm = cnorm + &
150: & PetscConj(solution(i,j)-userx(i,j))*(solution(i,j)-userx(i,j))
151: 80 continue
152: 90 continue
153: enorm = PetscRealPart(cnorm*hx*hy)
154: write(6,115) m,n,enorm
155: 115 format ('m = ',I2,' n = ',I2,' error norm = ',1PE11.4)
156: 100 continue
158: ! We are finished solving linear systems, so we clean up the
159: ! data structures.
161: DEALLOCATE (userx,userb,solution,rho)
163: call UserFinalizeLinearSolver(userctx,ierr)
164: call PetscFinalize(ierr)
165: end
167: ! ----------------------------------------------------------------
168: subroutine UserInitializeLinearSolver(m,n,userctx,ierr)
170: implicit none
172: #include <finclude/petscsys.h>
173: #include <finclude/petscvec.h>
174: #include <finclude/petscmat.h>
175: #include <finclude/petscksp.h>
176: #include <finclude/petscpc.h>
178: PetscInt m,n
179: PetscErrorCode ierr
180: PetscFortranAddr userctx(*)
182: common /param/ hx2,hy2
183: PetscReal hx2,hy2
185: ! Local variable declararions
186: Mat A
187: Vec b,x
188: KSP ksp
189: PetscInt Ntot,five,one
192: ! Here we assume use of a grid of size m x n, with all points on the
193: ! interior of the domain, i.e., we do not include the points corresponding
194: ! to homogeneous Dirichlet boundary conditions. We assume that the domain
195: ! is [0,1]x[0,1].
197: hx2 = (m+1)*(m+1)
198: hy2 = (n+1)*(n+1)
199: Ntot = m*n
201: five = 5
202: one = 1
204: ! Create the sparse matrix. Preallocate 5 nonzeros per row.
206: call MatCreateSeqAIJ(PETSC_COMM_SELF,Ntot,Ntot,five, &
207: & PETSC_NULL_INTEGER,A,ierr)
208: !
209: ! Create vectors. Here we create vectors with no memory allocated.
210: ! This way, we can use the data structures already in the program
211: ! by using VecPlaceArray() subroutine at a later stage.
212: !
213: call VecCreateSeqWithArray(PETSC_COMM_SELF,one,Ntot, &
214: & PETSC_NULL_SCALAR,b,ierr)
215: call VecDuplicate(b,x,ierr)
217: ! Create linear solver context. This will be used repeatedly for all
218: ! the linear solves needed.
220: call KSPCreate(PETSC_COMM_SELF,ksp,ierr)
222: userctx(1) = x
223: userctx(2) = b
224: userctx(3) = A
225: userctx(4) = ksp
226: userctx(5) = m
227: userctx(6) = n
229: return
230: end
231: ! -----------------------------------------------------------------------
233: ! Solves -div (rho grad psi) = F using finite differences.
234: ! rho is a 2-dimensional array of size m by n, stored in Fortran
235: ! style by columns. userb is a standard one-dimensional array.
237: subroutine UserDoLinearSolver(rho,userctx,userb,userx,ierr)
239: implicit none
241: #include <finclude/petscsys.h>
242: #include <finclude/petscvec.h>
243: #include <finclude/petscmat.h>
244: #include <finclude/petscksp.h>
245: #include <finclude/petscpc.h>
247: PetscErrorCode ierr
248: PetscFortranAddr userctx(*)
249: PetscScalar rho(*),userb(*),userx(*)
252: common /param/ hx2,hy2
253: PetscReal hx2,hy2
255: PC pc
256: KSP ksp
257: Vec b,x
258: Mat A
259: PetscInt m,n,one
260: PetscInt i,j,II,JJ
261: PetscScalar v
263: ! PetscScalar tmpx(*),tmpb(*)
265: one = 1
266: x = userctx(1)
267: b = userctx(2)
268: A = userctx(3)
269: ksp = userctx(4)
270: m = int(userctx(5))
271: n = int(userctx(6))
273: ! This is not the most efficient way of generating the matrix,
274: ! but let's not worry about it. We should have separate code for
275: ! the four corners, each edge and then the interior. Then we won't
276: ! have the slow if-tests inside the loop.
277: !
278: ! Compute the operator
279: ! -div rho grad
280: ! on an m by n grid with zero Dirichlet boundary conditions. The rho
281: ! is assumed to be given on the same grid as the finite difference
282: ! stencil is applied. For a staggered grid, one would have to change
283: ! things slightly.
285: II = 0
286: do 110 j=1,n
287: do 100 i=1,m
288: if (j .gt. 1) then
289: JJ = II - m
290: v = -0.5*(rho(II+1) + rho(JJ+1))*hy2
291: call MatSetValues(A,one,II,one,JJ,v,INSERT_VALUES,ierr)
292: endif
293: if (j .lt. n) then
294: JJ = II + m
295: v = -0.5*(rho(II+1) + rho(JJ+1))*hy2
296: call MatSetValues(A,one,II,one,JJ,v,INSERT_VALUES,ierr)
297: endif
298: if (i .gt. 1) then
299: JJ = II - 1
300: v = -0.5*(rho(II+1) + rho(JJ+1))*hx2
301: call MatSetValues(A,one,II,one,JJ,v,INSERT_VALUES,ierr)
302: endif
303: if (i .lt. m) then
304: JJ = II + 1
305: v = -0.5*(rho(II+1) + rho(JJ+1))*hx2
306: call MatSetValues(A,one,II,one,JJ,v,INSERT_VALUES,ierr)
307: endif
308: v = 2*rho(II+1)*(hx2+hy2)
309: call MatSetValues(A,one,II,one,II,v,INSERT_VALUES,ierr)
310: II = II+1
311: 100 continue
312: 110 continue
313: !
314: ! Assemble matrix
315: !
316: call MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY,ierr)
317: call MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY,ierr)
319: !
320: ! Set operators. Here the matrix that defines the linear system
321: ! also serves as the preconditioning matrix. Since all the matrices
322: ! will have the same nonzero pattern here, we indicate this so the
323: ! linear solvers can take advantage of this.
324: !
325: call KSPSetOperators(ksp,A,A,ierr)
327: !
328: ! Set linear solver defaults for this problem (optional).
329: ! - Here we set it to use direct LU factorization for the solution
330: !
331: call KSPGetPC(ksp,pc,ierr)
332: call PCSetType(pc,PCLU,ierr)
334: !
335: ! Set runtime options, e.g.,
336: ! -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
337: ! These options will override those specified above as long as
338: ! KSPSetFromOptions() is called _after_ any other customization
339: ! routines.
340: !
341: ! Run the program with the option -help to see all the possible
342: ! linear solver options.
343: !
344: call KSPSetFromOptions(ksp,ierr)
346: !
347: ! This allows the PETSc linear solvers to compute the solution
348: ! directly in the user's array rather than in the PETSc vector.
349: !
350: ! This is essentially a hack and not highly recommend unless you
351: ! are quite comfortable with using PETSc. In general, users should
352: ! write their entire application using PETSc vectors rather than
353: ! arrays.
354: !
355: call VecPlaceArray(x,userx,ierr)
356: call VecPlaceArray(b,userb,ierr)
358: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
359: ! Solve the linear system
360: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
362: call KSPSolve(ksp,b,x,ierr)
364: call VecResetArray(x,ierr)
365: call VecResetArray(b,ierr)
367: return
368: end
370: ! ------------------------------------------------------------------------
372: subroutine UserFinalizeLinearSolver(userctx,ierr)
374: implicit none
376: #include <finclude/petscsys.h>
377: #include <finclude/petscvec.h>
378: #include <finclude/petscmat.h>
379: #include <finclude/petscksp.h>
380: #include <finclude/petscpc.h>
382: PetscErrorCode ierr
383: PetscFortranAddr userctx(*)
385: ! Local variable declararions
387: Vec x,b
388: Mat A
389: KSP ksp
390: !
391: ! We are all done and don't need to solve any more linear systems, so
392: ! we free the work space. All PETSc objects should be destroyed when
393: ! they are no longer needed.
394: !
395: x = userctx(1)
396: b = userctx(2)
397: A = userctx(3)
398: ksp = userctx(4)
400: call VecDestroy(x,ierr)
401: call VecDestroy(b,ierr)
402: call MatDestroy(A,ierr)
403: call KSPDestroy(ksp,ierr)
405: return
406: end