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: ! petscsys.h - base PETSc routines
33: ! petscvec.h - vectors
34: ! petscmat.h - matrices
35: ! petscpc.h - preconditioners
36: ! petscksp.h - Krylov subspace methods
37: ! Additional include statements may be needed if using other PETSc
38: ! routines in a Fortran program, e.g.,
39: ! petscviewer.h - viewers
40: ! petscis.h - index sets
41: !
42: #include finclude/petscsys.h
43: #include finclude/petscvec.h
44: #include finclude/petscmat.h
45: #include finclude/petscpc.h
46: #include finclude/petscksp.h
47: !
48: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
49: ! Variable declarations
50: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
51: !
52: ! Variables:
53: ! ksp - linear solver context
54: ! x, b, u - approx solution, right-hand-side, exact solution vectors
55: ! A - matrix that defines linear system
56: ! its - iterations for convergence
57: ! norm - norm of error in solution
58: ! rctx - random number context
59: !
61: KSP ksp
62: Mat A
63: Vec x,b,u
64: PetscRandom rctx
65: double precision norm,h2,sigma1
66: PetscScalar none,sigma2,v,pfive,czero
67: PetscScalar cone
68: PetscInt dim,its,n,Istart
69: PetscInt Iend,i,j,II,JJ,one
70: PetscErrorCode ierr
71: PetscMPIInt rank
72: PetscTruth flg
73: logical use_random
75: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
76: ! Beginning of program
77: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
79: call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
80: #if !defined(PETSC_USE_COMPLEX)
81: write(6,*) "This example requires complex numbers."
82: goto 200
83: #endif
85: none = -1.0
86: n = 6
87: sigma1 = 100.0
88: czero = 0.0
89: cone = PETSC_i
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) 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: call PetscRandomSetInterval(rctx,czero,cone,ierr)
133: endif
134: h2 = 1.0/((n+1)*(n+1))
136: one = 1
137: do 10, II=Istart,Iend-1
138: v = -1.0
139: i = II/n
140: j = II - i*n
141: if (i.gt.0) then
142: JJ = II - n
143: call MatSetValues(A,one,II,one,JJ,v,ADD_VALUES,ierr)
144: endif
145: if (i.lt.n-1) then
146: JJ = II + n
147: call MatSetValues(A,one,II,one,JJ,v,ADD_VALUES,ierr)
148: endif
149: if (j.gt.0) then
150: JJ = II - 1
151: call MatSetValues(A,one,II,one,JJ,v,ADD_VALUES,ierr)
152: endif
153: if (j.lt.n-1) then
154: JJ = II + 1
155: call MatSetValues(A,one,II,one,JJ,v,ADD_VALUES,ierr)
156: endif
157: if (use_random) call PetscRandomGetValue(rctx, &
158: & sigma2,ierr)
159: v = 4.0 - sigma1*h2 + sigma2*h2
160: call MatSetValues(A,one,II,one,II,v,ADD_VALUES,ierr)
161: 10 continue
162: if (use_random) call PetscRandomDestroy(rctx,ierr)
164: ! Assemble matrix, using the 2-step process:
165: ! MatAssemblyBegin(), MatAssemblyEnd()
166: ! Computations can be done while messages are in transition
167: ! by placing code between these two statements.
169: call MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY,ierr)
170: call MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY,ierr)
172: ! Create parallel vectors.
173: ! - Here, the parallel partitioning of the vector is determined by
174: ! PETSc at runtime. We could also specify the local dimensions
175: ! if desired.
176: ! - Note: We form 1 vector from scratch and then duplicate as needed.
178: call VecCreate(PETSC_COMM_WORLD,u,ierr)
179: call VecSetSizes(u,PETSC_DECIDE,dim,ierr)
180: call VecSetFromOptions(u,ierr)
181: call VecDuplicate(u,b,ierr)
182: call VecDuplicate(b,x,ierr)
184: ! Set exact solution; then compute right-hand-side vector.
186: if (use_random) then
187: call PetscRandomCreate(PETSC_COMM_WORLD,rctx,ierr)
188: call PetscRandomSetFromOptions(rctx,ierr)
189: call VecSetRandom(u,rctx,ierr)
190: else
191: pfive = 0.5
192: call VecSet(u,pfive,ierr)
193: endif
194: call MatMult(A,u,b,ierr)
196: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
197: ! Create the linear solver and set various options
198: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
200: ! Create linear solver context
202: call KSPCreate(PETSC_COMM_WORLD,ksp,ierr)
204: ! Set operators. Here the matrix that defines the linear system
205: ! also serves as the preconditioning matrix.
207: call KSPSetOperators(ksp,A,A,DIFFERENT_NONZERO_PATTERN,ierr)
209: ! Set runtime options, e.g.,
210: ! -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
212: call KSPSetFromOptions(ksp,ierr)
214: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
215: ! Solve the linear system
216: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
218: call KSPSolve(ksp,b,x,ierr)
220: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
221: ! Check solution and clean up
222: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
224: ! Check the error
226: call VecAXPY(x,none,u,ierr)
227: call VecNorm(x,NORM_2,norm,ierr)
228: call KSPGetIterationNumber(ksp,its,ierr)
229: if (rank .eq. 0) then
230: if (norm .gt. 1.e-12) then
231: write(6,100) norm,its
232: else
233: write(6,110) its
234: endif
235: endif
236: 100 format('Norm of error ',e10.4,',iterations ',i5)
237: 110 format('Norm of error < 1.e-12,iterations ',i5)
239: ! Free work space. All PETSc objects should be destroyed when they
240: ! are no longer needed.
242: if (use_random) call PetscRandomDestroy(rctx,ierr)
243: call KSPDestroy(ksp,ierr)
244: call VecDestroy(u,ierr)
245: call VecDestroy(x,ierr)
246: call VecDestroy(b,ierr)
247: call MatDestroy(A,ierr)
249: #if !defined(PETSC_USE_COMPLEX)
250: 200 continue
251: #endif
252: call PetscFinalize(ierr)
253: end