Actual source code: ex11f.F
1: !
2: ! Description: Solves a complex linear system in parallel with KSP (Fortran code).
3: !
4: !/*T
5: ! Concepts: KSP^solving a Helmholtz equation
6: ! Concepts: complex numbers
7: ! Processors: n
8: !T*/
9: !
10: ! The model problem:
11: ! Solve Helmholtz equation on the unit square: (0,1) x (0,1)
12: ! -delta u - sigma1*u + i*sigma2*u = f,
13: ! where delta = Laplace operator
14: ! Dirichlet b.c.'s on all sides
15: ! Use the 2-D, five-point finite difference stencil.
16: !
17: ! Compiling the code:
18: ! This code uses the complex numbers version of PETSc, so configure
19: ! must be run to enable this
20: !
21: !
22: ! -----------------------------------------------------------------------
24: program main
25: implicit none
27: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
28: ! Include files
29: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
30: !
31: ! The following include statements are required for KSP Fortran programs:
32: ! petsc.h - base PETSc routines
33: ! petscvec.h - vectors
34: ! petscmat.h - matrices
35: ! petscpc.h - preconditioners
36: ! petscksp.h - Krylov subspace methods
37: ! Include the following to use PETSc random numbers:
38: ! petscsys.h - system routines
39: ! Additional include statements may be needed if using other PETSc
40: ! routines in a Fortran program, e.g.,
41: ! petscviewer.h - viewers
42: ! petscis.h - index sets
43: !
44: #include include/finclude/petsc.h
45: #include include/finclude/petscvec.h
46: #include include/finclude/petscmat.h
47: #include include/finclude/petscpc.h
48: #include include/finclude/petscksp.h
49: #include include/finclude/petscsys.h
50: !
51: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
52: ! Variable declarations
53: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
54: !
55: ! Variables:
56: ! ksp - linear solver context
57: ! x, b, u - approx solution, right-hand-side, exact solution vectors
58: ! A - matrix that defines linear system
59: ! its - iterations for convergence
60: ! norm - norm of error in solution
61: ! rctx - random number context
62: !
64: KSP ksp
65: Mat A
66: Vec x,b,u
67: PetscRandom rctx
68: double precision norm,h2,sigma1
69: PetscScalar none,sigma2,v,pfive
70: PetscInt dim,its,n,Istart
71: PetscInt Iend,i,j,II,JJ,one
72: PetscErrorCode ierr
73: PetscMPIInt rank
74: PetscTruth flg
75: logical use_random
77: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
78: ! Beginning of program
79: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
81: call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
82: #if !defined(PETSC_USE_COMPLEX)
83: write(6,*) "This example requires complex numbers."
84: goto 200
85: #endif
87: none = -1.0
88: n = 6
89: sigma1 = 100.0
90: call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
91: call PetscOptionsGetReal(PETSC_NULL_CHARACTER,'-sigma1',sigma1, &
92: & flg,ierr)
93: call PetscOptionsGetInt(PETSC_NULL_CHARACTER,'-n',n,flg,ierr)
94: dim = n*n
96: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
97: ! Compute the matrix and right-hand-side vector that define
98: ! the linear system, Ax = b.
99: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
101: ! Create parallel matrix, specifying only its global dimensions.
102: ! When using MatCreate(), the matrix format can be specified at
103: ! runtime. Also, the parallel partitioning of the matrix is
104: ! determined by PETSc at runtime.
106: call MatCreate(PETSC_COMM_WORLD,A,ierr)
107: call MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,dim,dim,ierr)
108: call MatSetFromOptions(A,ierr)
110: ! Currently, all PETSc parallel matrix formats are partitioned by
111: ! contiguous chunks of rows across the processors. Determine which
112: ! rows of the matrix are locally owned.
114: call MatGetOwnershipRange(A,Istart,Iend,ierr)
116: ! Set matrix elements in parallel.
117: ! - Each processor needs to insert only elements that it owns
118: ! locally (but any non-local elements will be sent to the
119: ! appropriate processor during matrix assembly).
120: ! - Always specify global rows and columns of matrix entries.
122: call PetscOptionsHasName(PETSC_NULL_CHARACTER,'-norandom', &
123: & flg,ierr)
124: if (flg .eq. 1) then
125: use_random = .false.
126: sigma2 = 10.0*PETSC_i
127: else
128: use_random = .true.
129: call PetscRandomCreate(PETSC_COMM_WORLD, &
130: & rctx,ierr)
131: call PetscRandomSetFromOptions(rctx,ierr)
132: endif
133: h2 = 1.0/((n+1)*(n+1))
135: one = 1
136: do 10, II=Istart,Iend-1
137: v = -1.0
138: i = II/n
139: j = II - i*n
140: if (i.gt.0) then
141: JJ = II - n
142: call MatSetValues(A,one,II,one,JJ,v,ADD_VALUES,ierr)
143: endif
144: if (i.lt.n-1) then
145: JJ = II + n
146: call MatSetValues(A,one,II,one,JJ,v,ADD_VALUES,ierr)
147: endif
148: if (j.gt.0) then
149: JJ = II - 1
150: call MatSetValues(A,one,II,one,JJ,v,ADD_VALUES,ierr)
151: endif
152: if (j.lt.n-1) then
153: JJ = II + 1
154: call MatSetValues(A,one,II,one,JJ,v,ADD_VALUES,ierr)
155: endif
156: if (use_random) call PetscRandomGetValueImaginary(rctx, &
157: & sigma2,ierr)
158: v = 4.0 - sigma1*h2 + sigma2*h2
159: call MatSetValues(A,one,II,one,II,v,ADD_VALUES,ierr)
160: 10 continue
161: if (use_random) call PetscRandomDestroy(rctx,ierr)
163: ! Assemble matrix, using the 2-step process:
164: ! MatAssemblyBegin(), MatAssemblyEnd()
165: ! Computations can be done while messages are in transition
166: ! by placing code between these two statements.
168: call MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY,ierr)
169: call MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY,ierr)
171: ! Create parallel vectors.
172: ! - Here, the parallel partitioning of the vector is determined by
173: ! PETSc at runtime. We could also specify the local dimensions
174: ! if desired.
175: ! - Note: We form 1 vector from scratch and then duplicate as needed.
177: call VecCreate(PETSC_COMM_WORLD,u,ierr)
178: call VecSetSizes(u,PETSC_DECIDE,dim,ierr)
179: call VecSetFromOptions(u,ierr)
180: call VecDuplicate(u,b,ierr)
181: call VecDuplicate(b,x,ierr)
183: ! Set exact solution; then compute right-hand-side vector.
185: if (use_random) then
186: call PetscRandomCreate(PETSC_COMM_WORLD,rctx,ierr)
187: call PetscRandomSetFromOptions(rctx,ierr)
188: call VecSetRandom(u,rctx,ierr)
189: else
190: pfive = 0.5
191: call VecSet(u,pfive,ierr)
192: endif
193: call MatMult(A,u,b,ierr)
195: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
196: ! Create the linear solver and set various options
197: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
199: ! Create linear solver context
201: call KSPCreate(PETSC_COMM_WORLD,ksp,ierr)
203: ! Set operators. Here the matrix that defines the linear system
204: ! also serves as the preconditioning matrix.
206: call KSPSetOperators(ksp,A,A,DIFFERENT_NONZERO_PATTERN,ierr)
208: ! Set runtime options, e.g.,
209: ! -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
211: call KSPSetFromOptions(ksp,ierr)
213: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
214: ! Solve the linear system
215: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
217: call KSPSolve(ksp,b,x,ierr)
219: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
220: ! Check solution and clean up
221: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
223: ! Check the error
225: call VecAXPY(x,none,u,ierr)
226: call VecNorm(x,NORM_2,norm,ierr)
227: call KSPGetIterationNumber(ksp,its,ierr)
228: if (rank .eq. 0) then
229: if (norm .gt. 1.e-12) then
230: write(6,100) norm,its
231: else
232: write(6,110) its
233: endif
234: endif
235: 100 format('Norm of error ',e10.4,',iterations ',i5)
236: 110 format('Norm of error < 1.e-12,iterations ',i5)
238: ! Free work space. All PETSc objects should be destroyed when they
239: ! are no longer needed.
241: if (use_random) call PetscRandomDestroy(rctx,ierr)
242: call KSPDestroy(ksp,ierr)
243: call VecDestroy(u,ierr)
244: call VecDestroy(x,ierr)
245: call VecDestroy(b,ierr)
246: call MatDestroy(A,ierr)
248: 200 continue
249: call PetscFinalize(ierr)
250: end