Actual source code: ex11f.F
1: !
2: ! "$Id: ex11f.F,v 1.31 2001/08/07 03:04:00 balay Exp $";
3: !
4: ! Description: Solves a complex linear system in parallel with SLES (Fortran code).
5: !
6: !/*T
7: ! Concepts: SLES^solving a Helmholtz equation
8: ! Concepts: complex numbers
9: ! Processors: n
10: !T*/
11: !
12: ! The model problem:
13: ! Solve Helmholtz equation on the unit square: (0,1) x (0,1)
14: ! -delta u - sigma1*u + i*sigma2*u = f,
15: ! where delta = Laplace operator
16: ! Dirichlet b.c.'s on all sides
17: ! Use the 2-D, five-point finite difference stencil.
18: !
19: ! Compiling the code:
20: ! This code uses the complex numbers version of PETSc, so one of the
21: ! following values of BOPT must be used for compiling the PETSc libraries
22: ! and this example:
23: ! BOPT=g_complex - debugging version
24: ! BOPT=O_complex - optimized version
25: ! BOPT=Opg_complex - profiling version
26: !
27: ! -----------------------------------------------------------------------
29: program main
30: implicit none
32: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
33: ! Include files
34: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
35: !
36: ! The following include statements are required for SLES Fortran programs:
37: ! petsc.h - base PETSc routines
38: ! petscvec.h - vectors
39: ! petscmat.h - matrices
40: ! petscpc.h - preconditioners
41: ! petscksp.h - Krylov subspace methods
42: ! petscsles.h - SLES interface
43: ! Include the following to use PETSc random numbers:
44: ! petscsys.h - system routines
45: ! Additional include statements may be needed if using other PETSc
46: ! routines in a Fortran program, e.g.,
47: ! petscviewer.h - viewers
48: ! petscis.h - index sets
49: !
50: #include include/finclude/petsc.h
51: #include include/finclude/petscvec.h
52: #include include/finclude/petscmat.h
53: #include include/finclude/petscpc.h
54: #include include/finclude/petscksp.h
55: #include include/finclude/petscsles.h
56: #include include/finclude/petscsys.h
57: !
58: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
59: ! Variable declarations
60: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
61: !
62: ! Variables:
63: ! sles - linear solver context
64: ! x, b, u - approx solution, right-hand-side, exact solution vectors
65: ! A - matrix that defines linear system
66: ! its - iterations for convergence
67: ! norm - norm of error in solution
68: ! rctx - random number context
69: !
71: SLES sles
72: Mat A
73: Vec x,b,u
74: PetscRandom rctx
75: double precision norm,h2,sigma1
76: PetscScalar none,sigma2,v,pfive
77: integer dim,flg,its,ierr,n,rank
78: integer Istart,Iend,i,j,II,JJ
79: logical use_random
81: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
82: ! Beginning of program
83: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
85: call PetscInitialize(PETSC_NULL_CHARACTER,ierr)
86: #if !defined(PETSC_USE_COMPLEX)
87: write(6,*) "This example requires complex numbers."
88: goto 200
89: #endif
91: none = -1.0
92: n = 6
93: sigma1 = 100.0
94: call MPI_Comm_rank(PETSC_COMM_WORLD,rank,ierr)
95: call PetscOptionsGetReal(PETSC_NULL_CHARACTER,'-sigma1',sigma1, &
96: & flg,ierr)
97: call PetscOptionsGetInt(PETSC_NULL_CHARACTER,'-n',n,flg,ierr)
98: dim = n*n
100: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
101: ! Compute the matrix and right-hand-side vector that define
102: ! the linear system, Ax = b.
103: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
105: ! Create parallel matrix, specifying only its global dimensions.
106: ! When using MatCreate(), the matrix format can be specified at
107: ! runtime. Also, the parallel partitioning of the matrix is
108: ! determined by PETSc at runtime.
110: call MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,dim, &
111: & dim,A,ierr)
112: call MatSetFromOptions(A,ierr)
114: ! Currently, all PETSc parallel matrix formats are partitioned by
115: ! contiguous chunks of rows across the processors. Determine which
116: ! rows of the matrix are locally owned.
118: call MatGetOwnershipRange(A,Istart,Iend,ierr)
120: ! Set matrix elements in parallel.
121: ! - Each processor needs to insert only elements that it owns
122: ! locally (but any non-local elements will be sent to the
123: ! appropriate processor during matrix assembly).
124: ! - Always specify global rows and columns of matrix entries.
126: call PetscOptionsHasName(PETSC_NULL_CHARACTER,'-norandom', &
127: & flg,ierr)
128: if (flg .eq. 1) then
129: use_random = .false.
130: sigma2 = 10.0*PETSC_i
131: else
132: use_random = .true.
133: call PetscRandomCreate(PETSC_COMM_WORLD, &
134: & RANDOM_DEFAULT_IMAGINARY,rctx,ierr)
135: endif
136: h2 = 1.0/((n+1)*(n+1))
138: do 10, II=Istart,Iend-1
139: v = -1.0
140: i = II/n
141: j = II - i*n
142: if (i.gt.0) then
143: JJ = II - n
144: call MatSetValues(A,1,II,1,JJ,v,ADD_VALUES,ierr)
145: endif
146: if (i.lt.n-1) then
147: JJ = II + n
148: call MatSetValues(A,1,II,1,JJ,v,ADD_VALUES,ierr)
149: endif
150: if (j.gt.0) then
151: JJ = II - 1
152: call MatSetValues(A,1,II,1,JJ,v,ADD_VALUES,ierr)
153: endif
154: if (j.lt.n-1) then
155: JJ = II + 1
156: call MatSetValues(A,1,II,1,JJ,v,ADD_VALUES,ierr)
157: endif
158: if (use_random) call PetscRandomGetValue(rctx,sigma2,ierr)
159: v = 4.0 - sigma1*h2 + sigma2*h2
160: call MatSetValues(A,1,II,1,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,RANDOM_DEFAULT, &
188: & rctx,ierr)
189: call VecSetRandom(rctx,u,ierr)
190: else
191: pfive = 0.5
192: call VecSet(pfive,u,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 SLESCreate(PETSC_COMM_WORLD,sles,ierr)
204: ! Set operators. Here the matrix that defines the linear system
205: ! also serves as the preconditioning matrix.
207: call SLESSetOperators(sles,A,A,DIFFERENT_NONZERO_PATTERN, &
208: & ierr)
210: ! Set runtime options, e.g.,
211: ! -ksp_type <type> -pc_type <type> -ksp_monitor -ksp_rtol <rtol>
213: call SLESSetFromOptions(sles,ierr)
215: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
216: ! Solve the linear system
217: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
219: call SLESSolve(sles,b,x,its,ierr)
221: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
222: ! Check solution and clean up
223: ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
225: ! Check the error
227: call VecAXPY(none,u,x,ierr)
228: call VecNorm(x,NORM_2,norm,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 SLESDestroy(sles,ierr)
244: call VecDestroy(u,ierr)
245: call VecDestroy(x,ierr)
246: call VecDestroy(b,ierr)
247: call MatDestroy(A,ierr)
249: 200 continue
250: call PetscFinalize(ierr)
251: end