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