download the original source code.
1 c
2 c Example 5
3 c
4 c Interface: Linear-Algebraic (IJ), Fortran (77) version
5 c
6 c Compile with: make ex5f
7 c
8 c Sample run: mpirun -np 4 ex5f
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,
18 c and PCG with ParaSails
19 c
20 c
21 c Notes: for PCG, GMRES and BiCGStab, precond_id means:
22 c 0 - do not set up a preconditioner
23 c 1 - set up a ds preconditioner
24 c 2 - set up an amg preconditioner
25 c 3 - set up a pilut preconditioner
26 c 4 - set up a ParaSails preconditioner
27 c
28
29 program ex5f
30
31
32 implicit none
33
34 include 'mpif.h'
35
36 integer MAX_LOCAL_SIZE
37 integer HYPRE_PARCSR
38
39 parameter (MAX_LOCAL_SIZE=123000)
40
41 c the following is from HYPRE.c
42 parameter (HYPRE_PARCSR=5555)
43
44 integer ierr
45 integer num_procs, myid
46 integer local_size, extra
47 integer n, solver_id, print_solution, ng
48 integer nnz, ilower, iupper, i
49 integer precond_id;
50 double precision h, h2
51 double precision rhs_values(MAX_LOCAL_SIZE)
52 double precision x_values(MAX_LOCAL_SIZE)
53 integer rows(MAX_LOCAL_SIZE)
54 integer cols(5)
55 double precision values(5)
56 integer num_iterations
57 double precision final_res_norm, tol
58
59 integer*8 mpi_comm
60 integer*8 parcsr_A
61 integer*8 A
62 integer*8 b
63 integer*8 x
64 integer*8 par_b
65 integer*8 par_x
66 integer*8 solver
67 integer*8 precond
68
69 c-----------------------------------------------------------------------
70 c Initialize MPI
71 c-----------------------------------------------------------------------
72
73 call MPI_INIT(ierr)
74 call MPI_COMM_RANK(MPI_COMM_WORLD, myid, ierr)
75 call MPI_COMM_SIZE(MPI_COMM_WORLD, num_procs, ierr)
76 mpi_comm = MPI_COMM_WORLD
77
78 c Default problem parameters
79 n = 33
80 solver_id = 0
81 print_solution = 0
82 tol = 1.0d-7
83
84 c The input section not implemented yet.
85
86 c Preliminaries: want at least one processor per row
87 if ( n*n .lt. num_procs) then
88 n = int(sqrt(real(num_procs))) + 1
89 endif
90 c ng = global no. rows, h = mesh size
91 ng = n*n
92 h = 1.0d0/(n+1)
93 h2 = h*h
94
95 c Each processor knows only of its own rows - the range is denoted by ilower
96 c and upper. Here we partition the rows. We account for the fact that
97 c N may not divide evenly by the number of processors.
98 local_size = ng/num_procs
99 extra = ng - local_size*num_procs
100
101 ilower = local_size*myid
102 ilower = ilower + min(myid, extra)
103
104 iupper = local_size*(myid+1)
105 iupper = iupper + min(myid+1, extra)
106 iupper = iupper - 1
107
108 c How many rows do I have?
109 local_size = iupper - ilower + 1
110
111 c Create the matrix.
112 c Note that this is a square matrix, so we indicate the row partition
113 c size twice (since number of rows = number of cols)
114 call HYPRE_IJMatrixCreate( mpi_comm, ilower,
115 1 iupper, ilower, iupper, A, ierr )
116
117
118 c Choose a parallel csr format storage (see the User's Manual)
119 call HYPRE_IJMatrixSetObjectType( A, HYPRE_PARCSR, ierr)
120
121 c Initialize before setting coefficients
122 call HYPRE_IJMatrixInitialize( A, ierr)
123
124
125 c Now go through my local rows and set the matrix entries.
126 c Each row has at most 5 entries. For example, if n=3:
127 c
128 c A = [M -I 0; -I M -I; 0 -I M]
129 c M = [4 -1 0; -1 4 -1; 0 -1 4]
130 c
131 c Note that here we are setting one row at a time, though
132 c one could set all the rows together (see the User's Manual).
133
134
135 do i = ilower, iupper
136 nnz = 1
137
138
139 c The left identity block:position i-n
140 if ( (i-n) .ge. 0 ) then
141 cols(nnz) = i-n
142 values(nnz) = -1.0d0
143 nnz = nnz + 1
144 endif
145
146 c The left -1: position i-1
147 if ( mod(i,n).ne.0 ) then
148 cols(nnz) = i-1
149 values(nnz) = -1.0d0
150 nnz = nnz + 1
151 endif
152
153 c Set the diagonal: position i
154 cols(nnz) = i
155 values(nnz) = 4.0d0
156 nnz = nnz + 1
157
158 c The right -1: position i+1
159 if ( mod((i+1),n) .ne. 0 ) then
160 cols(nnz) = i+1
161 values(nnz) = -1.0d0
162 nnz = nnz + 1
163 endif
164
165 c The right identity block:position i+n
166 if ((i+n) .lt. ng ) then
167 cols(nnz) = i+n
168 values(nnz) = -1.0d0
169 nnz = nnz + 1
170 endif
171
172 c Set the values for row i
173 call HYPRE_IJMatrixSetValues(
174 1 A, 1, nnz-1, i, cols, values, ierr)
175
176 enddo
177
178
179 c Assemble after setting the coefficients
180 call HYPRE_IJMatrixAssemble( A, ierr)
181
182 c Get parcsr matrix object
183 call HYPRE_IJMatrixGetObject( A, parcsr_A, ierr)
184
185
186 c Create the rhs and solution
187 call HYPRE_IJVectorCreate(mpi_comm,
188 1 ilower, iupper, b, ierr )
189 call HYPRE_IJVectorSetObjectType(b, HYPRE_PARCSR, ierr)
190 call HYPRE_IJVectorInitialize(b, ierr)
191
192 call HYPRE_IJVectorCreate(mpi_comm,
193 1 ilower, iupper, x, ierr )
194 call HYPRE_IJVectorSetObjectType(x, HYPRE_PARCSR, ierr)
195 call HYPRE_IJVectorInitialize(x, ierr)
196
197
198 c Set the rhs values to h^2 and the solution to zero
199 do i = 1, local_size
200 rhs_values(i) = h2
201 x_values(i) = 0.0
202 rows(i) = ilower + i -1
203 enddo
204 call HYPRE_IJVectorSetValues(
205 1 b, local_size, rows, rhs_values, ierr )
206 call HYPRE_IJVectorSetValues(
207 1 x, local_size, rows, x_values, ierr)
208
209
210 call HYPRE_IJVectorAssemble( b, ierr)
211 call HYPRE_IJVectorAssemble( x, ierr)
212
213 c get the x and b objects
214
215 call HYPRE_IJVectorGetObject( b, par_b, ierr)
216 call HYPRE_IJVectorGetObject( x, par_x, ierr)
217
218
219 c Choose a solver and solve the system
220
221 c AMG
222 if ( solver_id .eq. 0 ) then
223
224 c Create solver
225 call HYPRE_BoomerAMGCreate(solver, ierr)
226
227
228 c Set some parameters (See Reference Manual for more parameters)
229
230 c print solve info + parameters
231 call HYPRE_BoomerAMGSetPrintLevel(solver, 3, ierr)
232 c Falgout coarsening
233 call HYPRE_BoomerAMGSetCoarsenType(solver, 6, ierr)
234 c G-S/Jacobi hybrid relaxation
235 call HYPRE_BoomerAMGSetRelaxType(solver, 3, ierr)
236 c Sweeeps on each level
237 call HYPRE_BoomerAMGSetNumSweeps(solver, 1, ierr)
238 c maximum number of levels
239 call HYPRE_BoomerAMGSetMaxLevels(solver, 20, ierr)
240 c conv. tolerance
241 call HYPRE_BoomerAMGSetTol(solver, 1.0d-7, ierr)
242
243 c Now setup and solve!
244 call HYPRE_BoomerAMGSetup(
245 1 solver, parcsr_A, par_b, par_x, ierr )
246 call HYPRE_BoomerAMGSolve(
247 1 solver, parcsr_A, par_b, par_x, ierr )
248
249
250 c Run info - needed logging turned on
251 call HYPRE_BoomerAMGGetNumIterations(solver, num_iterations,
252 1 ierr)
253 call HYPRE_BoomerAMGGetFinalReltvRes(solver, final_res_norm,
254 1 ierr)
255
256
257 if (myid .eq. 0) then
258 print *
259 print '(A,I2)', " Iterations = ", num_iterations
260 print '(A,ES16.8)',
261 1 " Final Relative Residual Norm = ", final_res_norm
262 print *
263 endif
264
265 c Destroy solver
266 call HYPRE_BoomerAMGDestroy( solver, ierr )
267
268 c PCG (with DS)
269 elseif (solver_id .eq. 50) then
270
271
272 c Create solver
273 call HYPRE_ParCSRPCGCreate(MPI_COMM_WORLD, solver, ierr)
274
275 c Set some parameters (See Reference Manual for more parameters)
276 call HYPRE_ParCSRPCGSetMaxIter(solver, 1000, ierr)
277 call HYPRE_ParCSRPCGSetTol(solver, 1.0d-7, ierr)
278 call HYPRE_ParCSRPCGSetTwoNorm(solver, 1, ierr)
279 call HYPRE_ParCSRPCGSetPrintLevel(solver, 2, ierr)
280 call HYPRE_ParCSRPCGSetLogging(solver, 1, ierr)
281
282 c set ds (diagonal scaling) as the pcg preconditioner
283 precond_id = 1
284 call HYPRE_ParCSRPCGSetPrecond(solver, precond_id,
285 1 precond, ierr)
286
287
288
289 c Now setup and solve!
290 call HYPRE_ParCSRPCGSetup(solver, parcsr_A, par_b,
291 & par_x, ierr)
292 call HYPRE_ParCSRPCGSolve(solver, parcsr_A, par_b,
293 & par_x, ierr)
294
295
296 c Run info - needed logging turned on
297
298 call HYPRE_ParCSRPCGGetNumIterations(solver, num_iterations,
299 & ierr)
300 call HYPRE_ParCSRPCGGetFinalRelative(solver, final_res_norm,
301 & ierr)
302
303 if (myid .eq. 0) then
304 print *
305 print *, "Iterations = ", num_iterations
306 print *, "Final Relative Residual Norm = ", final_res_norm
307 print *
308 endif
309
310 c Destroy solver
311 call HYPRE_ParCSRPCGDestroy(solver, ierr)
312
313
314 c PCG with AMG preconditioner
315 elseif (solver_id == 1) then
316
317 c Create solver
318 call HYPRE_ParCSRPCGCreate(MPI_COMM_WORLD, solver, ierr)
319
320 c Set some parameters (See Reference Manual for more parameters)
321 call HYPRE_ParCSRPCGSetMaxIter(solver, 1000, ierr)
322 call HYPRE_ParCSRPCGSetTol(solver, 1.0d-7, ierr)
323 call HYPRE_ParCSRPCGSetTwoNorm(solver, 1, ierr)
324 call HYPRE_ParCSRPCGSetPrintLevel(solver, 2, ierr)
325 call HYPRE_ParCSRPCGSetLogging(solver, 1, ierr)
326
327 c Now set up the AMG preconditioner and specify any parameters
328
329 call HYPRE_BoomerAMGCreate(precond, ierr)
330
331
332 c Set some parameters (See Reference Manual for more parameters)
333
334 c print less solver info since a preconditioner
335 call HYPRE_BoomerAMGSetPrintLevel(precond, 1, ierr);
336 c Falgout coarsening
337 call HYPRE_BoomerAMGSetCoarsenType(precond, 6, ierr)
338 c SYMMETRIC G-S/Jacobi hybrid relaxation
339 call HYPRE_BoomerAMGSetRelaxType(precond, 6, ierr)
340 c Sweeeps on each level
341 call HYPRE_BoomerAMGSetNumSweeps(precond, 1, ierr)
342 c conv. tolerance
343 call HYPRE_BoomerAMGSetTol(precond, 0.0d0, ierr)
344 c do only one iteration!
345 call HYPRE_BoomerAMGSetMaxIter(precond, 1, ierr)
346
347 c set amg as the pcg preconditioner
348 precond_id = 2
349 call HYPRE_ParCSRPCGSetPrecond(solver, precond_id,
350 1 precond, ierr)
351
352
353 c Now setup and solve!
354 call HYPRE_ParCSRPCGSetup(solver, parcsr_A, par_b,
355 1 par_x, ierr)
356 call HYPRE_ParCSRPCGSolve(solver, parcsr_A, par_b,
357 1 par_x, ierr)
358
359
360 c Run info - needed logging turned on
361
362 call HYPRE_ParCSRPCGGetNumIterations(solver, num_iterations,
363 1 ierr)
364 call HYPRE_ParCSRPCGGetFinalRelative(solver, final_res_norm,
365 1 ierr)
366
367 if (myid .eq. 0) then
368 print *
369 print *, "Iterations = ", num_iterations
370 print *, "Final Relative Residual Norm = ", final_res_norm
371 print *
372 endif
373
374 c Destroy precond and solver
375
376 call HYPRE_BoomerAMGDestroy(precond, ierr )
377 call HYPRE_ParCSRPCGDestroy(solver, ierr)
378
379 c PCG with ParaSails
380 elseif (solver_id .eq. 8) then
381
382 c Create solver
383 call HYPRE_ParCSRPCGCreate(MPI_COMM_WORLD, solver, ierr)
384
385 c Set some parameters (See Reference Manual for more parameters)
386 call HYPRE_ParCSRPCGSetMaxIter(solver, 1000, ierr)
387 call HYPRE_ParCSRPCGSetTol(solver, 1.0d-7, ierr)
388 call HYPRE_ParCSRPCGSetTwoNorm(solver, 1, ierr)
389 call HYPRE_ParCSRPCGSetPrintLevel(solver, 2, ierr)
390 call HYPRE_ParCSRPCGSetLogging(solver, 1, ierr)
391
392 c Now set up the Parasails preconditioner and specify any parameters
393 call HYPRE_ParaSailsCreate(MPI_COMM_WORLD, precond,ierr)
394 call HYPRE_ParaSailsSetParams(precond, 0.1d0, 1, ierr)
395 call HYPRE_ParaSailsSetFilter(precond, 0.05d0, ierr)
396 call HYPRE_ParaSailsSetSym(precond, 1)
397 call HYPRE_ParaSailsSetLogging(precond, 3, ierr)
398
399 c set parsails as the pcg preconditioner
400 precond_id = 4
401 call HYPRE_ParCSRPCGSetPrecond(solver, precond_id,
402 1 precond, ierr)
403
404
405 c Now setup and solve!
406 call HYPRE_ParCSRPCGSetup(solver, parcsr_A, par_b,
407 1 par_x, ierr)
408 call HYPRE_ParCSRPCGSolve(solver, parcsr_A, par_b,
409 1 par_x, ierr)
410
411
412 c Run info - needed logging turned on
413
414 call HYPRE_ParCSRPCGGetNumIterations(solver, num_iterations,
415 1 ierr)
416 call HYPRE_ParCSRPCGGetFinalRelative(solver, final_res_norm,
417 1 ierr)
418
419 if (myid .eq. 0) then
420 print *
421 print *, "Iterations = ", num_iterations
422 print *, "Final Relative Residual Norm = ", final_res_norm
423 print *
424 endif
425
426 c Destroy precond and solver
427
428 call HYPRE_ParaSailsDestroy(precond, ierr )
429 call HYPRE_ParCSRPCGDestroy(solver, ierr)
430
431 else
432 if (myid .eq. 0) then
433 print *,'Invalid solver id specified'
434 stop
435 endif
436 endif
437
438
439
440 c Print the solution
441 if ( print_solution .ne. 0 ) then
442 call HYPRE_IJVectorPrint( x, "ij.out.x", ierr)
443 endif
444
445 c Clean up
446
447 call HYPRE_IJMatrixDestroy(A, ierr)
448 call HYPRE_IJVectorDestroy(b, ierr)
449 call HYPRE_IJVectorDestroy(x, ierr)
450
451
452 c Finalize MPI
453 call MPI_Finalize(ierr)
454
455 stop
456 end
syntax highlighted by Code2HTML, v. 0.9.1