download the original source code.
1 c
2 c Example 5
3 c
4 c Interface: Linear-Algebraic (IJ), Babel-based version in Fortran
5 c
6 c Compile with: make ex5b
7 c
8 c Sample run: mpirun -np 4 ex5b
9 c
10 c Description: This example solves the 2-D
11 c Laplacian problem with zero boundary conditions
12 c on an nxn grid. The number of unknowns is N=n^2.
13 c The standard 5-point stencil is used, and we solve
14 c for the interior nodes only.
15 c
16 c This example solves the same problem as Example 3.
17 c Available solvers are AMG, PCG, and PCG with AMG or
18 c Parasails preconditioners.
19
20 program ex5b77
21
22
23 implicit none
24
25 include 'mpif.h'
26
27 integer MAX_LOCAL_SIZE
28 parameter (MAX_LOCAL_SIZE=123000)
29
30 integer ierr, ierrtmp
31 integer num_procs, myid
32 integer local_size, extra
33 integer n, solver_id, print_solution, ng
34 integer nnz, ilower, iupper, i
35 double precision h, h2
36 double precision rhs_values(MAX_LOCAL_SIZE)
37 double precision x_values(MAX_LOCAL_SIZE)
38 integer rows(MAX_LOCAL_SIZE)
39 integer cols(5)
40 double precision values(5)
41 integer num_iterations
42 double precision final_res_norm, tol
43
44 integer*8 bHYPRE_mpicomm
45 integer*8 mpi_comm
46 integer*8 parcsr_A
47 integer*8 op_A
48 integer*8 par_b
49 integer*8 par_x
50 integer*8 vec_b
51 integer*8 vec_x
52 integer*8 amg_solver
53 integer*8 except
54 c ... except is for Babel exceptions, which we shall ignore
55
56 c-----------------------------------------------------------------------
57 c Initialize MPI
58 c-----------------------------------------------------------------------
59
60 call MPI_INIT(ierr)
61 call MPI_COMM_RANK(MPI_COMM_WORLD, myid, ierr)
62 call MPI_COMM_SIZE(MPI_COMM_WORLD, num_procs, ierr)
63 mpi_comm = MPI_COMM_WORLD
64 call bHYPRE_MPICommunicator_CreateF_f( mpi_comm, bHYPRE_mpicomm,
65 1 except )
66
67 c Default problem parameters
68 n = 33
69 solver_id = 0
70 print_solution = 0
71 tol = 1.0d-7
72
73 c The input section not implemented yet.
74
75 c Preliminaries: want at least one processor per row
76 if ( n*n .lt. num_procs) then
77 n = int(sqrt(real(num_procs))) + 1
78 endif
79 c ng = global no. rows, h = mesh size
80 ng = n*n
81 h = 1.0d0/(n+1)
82 h2 = h*h
83
84 c Each processor knows only of its own rows - the range is denoted by ilower
85 c and upper. Here we partition the rows. We account for the fact that
86 c N may not divide evenly by the number of processors.
87 local_size = ng/num_procs
88 extra = ng - local_size*num_procs
89
90 ilower = local_size*myid
91 ilower = ilower + min(myid, extra)
92
93 iupper = local_size*(myid+1)
94 iupper = iupper + min(myid+1, extra)
95 iupper = iupper - 1
96
97 c How many rows do I have?
98 local_size = iupper - ilower + 1
99
100 c Create the matrix.
101 c Note that this is a square matrix, so we indicate the row partition
102 c size twice (since number of rows = number of cols)
103 call bHYPRE_IJParCSRMatrix_Create_f( bHYPRE_mpicomm, ilower,
104 1 iupper, ilower, iupper, parcsr_A, except )
105
106 c op_A will be needed later as a function argument
107 call bHYPRE_Operator__cast_f( parcsr_A, op_A, except )
108
109 c Choose a parallel csr format storage (see the User's Manual)
110 c Note: Here the HYPRE interface requires a SetObjectType call.
111 c I am using the bHYPRE interface in a way which does not because
112 c the object type is already specified through the class name.
113
114 c Initialize before setting coefficients
115 call bHYPRE_IJParCSRMatrix_Initialize_f( parcsr_A, ierrtmp,
116 1 except )
117
118 c Now go through my local rows and set the matrix entries.
119 c Each row has at most 5 entries. For example, if n=3:
120 c
121 c A = [M -I 0; -I M -I; 0 -I M]
122 c M = [4 -1 0; -1 4 -1; 0 -1 4]
123 c
124 c Note that here we are setting one row at a time, though
125 c one could set all the rows together (see the User's Manual).
126
127 do i = ilower, iupper
128 nnz = 1
129
130 c The left identity block:position i-n
131 if ( (i-n) .ge. 0 ) then
132 cols(nnz) = i-n
133 values(nnz) = -1.0d0
134 nnz = nnz + 1
135 endif
136
137 c The left -1: position i-1
138 if ( mod(i,n).ne.0 ) then
139 cols(nnz) = i-1
140 values(nnz) = -1.0d0
141 nnz = nnz + 1
142 endif
143
144 c Set the diagonal: position i
145 cols(nnz) = i
146 values(nnz) = 4.0d0
147 nnz = nnz + 1
148
149 c The right -1: position i+1
150 if ( mod((i+1),n) .ne. 0 ) then
151 cols(nnz) = i+1
152 values(nnz) = -1.0d0
153 nnz = nnz + 1
154 endif
155
156 c The right identity block:position i+n
157 if ((i+n) .lt. ng ) then
158 cols(nnz) = i+n
159 values(nnz) = -1.0d0
160 nnz = nnz + 1
161 endif
162
163 c Set the values for row i
164 call bHYPRE_IJParCSRMatrix_SetValues_f(
165 1 parcsr_A, 1, nnz-1, i, cols, values, 5, ierrtmp, except )
166
167 enddo
168
169
170 c Assemble after setting the coefficients
171 call bHYPRE_IJParCSRMatrix_Assemble_f( parcsr_A, ierrtmp, except )
172
173 c Create the rhs and solution
174 call bHYPRE_IJParCSRVector_Create_f( bHYPRE_mpicomm,
175 1 ilower, iupper, par_b, except )
176 c vec_b will be needed later for function arguments
177 call bHYPRE_Vector__cast_f( par_b, vec_b, except )
178
179 call bHYPRE_IJParCSRVector_Initialize_f( par_b, ierrtmp, except )
180
181 call bHYPRE_IJParCSRVector_Create_f( bHYPRE_mpicomm,
182 1 ilower, iupper, par_x, except )
183 c vec_x will be needed later for function arguments
184 call bHYPRE_Vector__cast_f( par_x, vec_x, except )
185
186 call bHYPRE_IJParCSRVector_Initialize_f( par_x, ierrtmp, except )
187
188 c Set the rhs values to h^2 and the solution to zero
189 do i = 1, local_size
190 rhs_values(i) = h2
191 x_values(i) = 0.0
192 rows(i) = ilower + i -1
193 enddo
194 call bHYPRE_IJParCSRVector_SetValues_f(
195 1 par_b, local_size, rows, rhs_values, ierrtmp, except )
196 call bHYPRE_IJParCSRVector_SetValues_f(
197 1 par_x, local_size, rows, x_values, ierrtmp, except )
198
199
200 call bHYPRE_IJParCSRVector_Assemble_f( par_b, ierrtmp, except )
201 call bHYPRE_IJParCSRVector_Assemble_f( par_x, ierrtmp, except )
202
203 c Choose a solver and solve the system
204
205 c AMG
206 if ( solver_id == 0 ) then
207
208 c Create solver
209 call bHYPRE_BoomerAMG_Create_f(
210 1 bHYPRE_mpicomm, parcsr_A, amg_solver, except )
211
212 c Set some parameters (See Reference Manual for more parameters)
213 c PrintLevel=3 means print solve info + parameters
214 c CoarsenType=6 means Falgout coarsening
215 c RelaxType=3 means Gauss-Seidel/Jacobi hybrid relaxation
216 call bHYPRE_BoomerAMG_SetIntParameter_f(
217 1 amg_solver, "PrintLevel", 3, ierrtmp, except )
218 call bHYPRE_BoomerAMG_SetIntParameter_f(
219 1 amg_solver, "CoarsenType", 6, ierrtmp, except )
220 call bHYPRE_BoomerAMG_SetIntParameter_f(
221 1 amg_solver, "RelaxType", 3, ierrtmp, except )
222 call bHYPRE_BoomerAMG_SetIntParameter_f(
223 1 amg_solver, "NumSweeps", 1, ierrtmp, except )
224 call bHYPRE_BoomerAMG_SetIntParameter_f(
225 1 amg_solver, "MaxLevels", 20, ierrtmp, except )
226 call bHYPRE_BoomerAMG_SetDoubleParameter_f(
227 1 amg_solver, "Tolerance", tol, ierrtmp, except )
228
229 c Now setup and solve!
230 call bHYPRE_BoomerAMG_Setup_f(
231 1 amg_solver, vec_b, vec_x, ierrtmp, except )
232 call bHYPRE_BoomerAMG_Apply_f(
233 1 amg_solver, vec_b, vec_x, ierrtmp, except )
234
235 c Run info - needed logging turned on
236 call bHYPRE_BoomerAMG_GetIntValue_f(
237 1 amg_solver, "NumIterations", num_iterations, ierrtmp,
238 2 except )
239 ierr = ierr + ierrtmp
240 call bHYPRE_BoomerAMG_GetDoubleValue_f(
241 1 amg_solver, "RelResidualNorm", final_res_norm, ierrtmp,
242 2 except )
243
244 if (myid .eq. 0) then
245 print *
246 print *, "Iterations = ", num_iterations
247 print *, "Final Relative Residual Norm = ", final_res_norm
248 print *
249 endif
250
251 c Destroy solver
252 call bHYPRE_BoomerAMG_deleteRef_f( amg_solver, except )
253
254 endif
255
256 c The calls of other solvers are not implemented yet.
257
258
259 c Print the solution
260 if ( print_solution .ne. 0 ) then
261 call bHYPRE_IJParCSRVector_Print_f( par_x, "ij.out.x", except )
262 endif
263
264 c Clean up
265 call bHYPRE_Operator_deleteRef_f( op_A, except )
266 call bHYPRE_Vector_deleteRef_f( vec_b, except )
267 call bHYPRE_Vector_deleteRef_f( vec_x, except )
268 call bHYPRE_IJParCSRMatrix_deleteRef_f( parcsr_A, except )
269 call bHYPRE_IJParCSRVector_deleteRef_f( par_b, except )
270 call bHYPRE_IJParCSRVector_deleteRef_f( par_x, except )
271 call bHYPRE_MPICommunicator_deleteRef_f( bHYPRE_mpicomm, except )
272
273 c We need a multi-language equivalent of hypre_assert.
274 if ( ierr .ne. 0 ) then
275 print *
276 print *, "***** Bad ierr = ", ierr
277 print *
278 endif
279
280 c Finalize MPI
281 call MPI_Finalize(ierrtmp)
282
283 stop
284 end
syntax highlighted by Code2HTML, v. 0.9.1