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 Falgout coarsening
234 call HYPRE_BoomerAMGSetCoarsenType(solver, 6, ierr)
235 c G-S/Jacobi hybrid relaxation
236 call HYPRE_BoomerAMGSetRelaxType(solver, 3, ierr)
237 c Sweeeps on each level
238 call HYPRE_BoomerAMGSetNumSweeps(solver, 1, ierr)
239 c maximum number of levels
240 call HYPRE_BoomerAMGSetMaxLevels(solver, 20, ierr)
241 c conv. tolerance
242 call HYPRE_BoomerAMGSetTol(solver, 1.0d-7, ierr)
243
244 c Now setup and solve!
245 call HYPRE_BoomerAMGSetup(
246 1 solver, parcsr_A, par_b, par_x, ierr)
247 call HYPRE_BoomerAMGSolve(
248 1 solver, parcsr_A, par_b, par_x, ierr)
249
250
251 c Run info - needed logging turned on
252 call HYPRE_BoomerAMGGetNumIterations(solver, num_iterations,
253 1 ierr)
254 call HYPRE_BoomerAMGGetFinalReltvRes(solver, final_res_norm,
255 1 ierr)
256
257
258 if ( myid .eq. 0 ) then
259 print *
260 print '(A,I2)', " Iterations = ", num_iterations
261 print '(A,ES16.8)',
262 1 " Final Relative Residual Norm = ", final_res_norm
263 print *
264 endif
265
266 c Destroy solver
267 call HYPRE_BoomerAMGDestroy(solver, ierr)
268
269 c PCG (with DS)
270 elseif ( solver_id .eq. 50 ) then
271
272
273 c Create solver
274 call HYPRE_ParCSRPCGCreate(MPI_COMM_WORLD, solver, ierr)
275
276 c Set some parameters (See Reference Manual for more parameters)
277 call HYPRE_ParCSRPCGSetMaxIter(solver, 1000, ierr)
278 call HYPRE_ParCSRPCGSetTol(solver, 1.0d-7, ierr)
279 call HYPRE_ParCSRPCGSetTwoNorm(solver, 1, ierr)
280 call HYPRE_ParCSRPCGSetPrintLevel(solver, 2, ierr)
281 call HYPRE_ParCSRPCGSetLogging(solver, 1, ierr)
282
283 c set ds (diagonal scaling) as the pcg preconditioner
284 precond_id = 1
285 call HYPRE_ParCSRPCGSetPrecond(solver, precond_id,
286 1 precond, ierr)
287
288
289
290 c Now setup and solve!
291 call HYPRE_ParCSRPCGSetup(solver, parcsr_A, par_b,
292 & par_x, ierr)
293 call HYPRE_ParCSRPCGSolve(solver, parcsr_A, par_b,
294 & par_x, ierr)
295
296
297 c Run info - needed logging turned on
298
299 call HYPRE_ParCSRPCGGetNumIterations(solver, num_iterations,
300 & ierr)
301 call HYPRE_ParCSRPCGGetFinalRelative(solver, final_res_norm,
302 & ierr)
303
304 if ( myid .eq. 0 ) then
305 print *
306 print *, "Iterations = ", num_iterations
307 print *, "Final Relative Residual Norm = ", final_res_norm
308 print *
309 endif
310
311 c Destroy solver
312 call HYPRE_ParCSRPCGDestroy(solver, ierr)
313
314
315 c PCG with AMG preconditioner
316 elseif ( solver_id == 1 ) then
317
318 c Create solver
319 call HYPRE_ParCSRPCGCreate(MPI_COMM_WORLD, solver, ierr)
320
321 c Set some parameters (See Reference Manual for more parameters)
322 call HYPRE_ParCSRPCGSetMaxIter(solver, 1000, ierr)
323 call HYPRE_ParCSRPCGSetTol(solver, 1.0d-7, ierr)
324 call HYPRE_ParCSRPCGSetTwoNorm(solver, 1, ierr)
325 call HYPRE_ParCSRPCGSetPrintLevel(solver, 2, ierr)
326 call HYPRE_ParCSRPCGSetLogging(solver, 1, ierr)
327
328 c Now set up the AMG preconditioner and specify any parameters
329
330 call HYPRE_BoomerAMGCreate(precond, ierr)
331
332
333 c Set some parameters (See Reference Manual for more parameters)
334
335 c print less solver info since a preconditioner
336 call HYPRE_BoomerAMGSetPrintLevel(precond, 1, ierr);
337 c Falgout coarsening
338 call HYPRE_BoomerAMGSetCoarsenType(precond, 6, ierr)
339 c SYMMETRIC G-S/Jacobi hybrid relaxation
340 call HYPRE_BoomerAMGSetRelaxType(precond, 6, ierr)
341 c Sweeeps on each level
342 call HYPRE_BoomerAMGSetNumSweeps(precond, 1, ierr)
343 c conv. tolerance
344 call HYPRE_BoomerAMGSetTol(precond, 0.0d0, ierr)
345 c do only one iteration!
346 call HYPRE_BoomerAMGSetMaxIter(precond, 1, ierr)
347
348 c set amg as the pcg preconditioner
349 precond_id = 2
350 call HYPRE_ParCSRPCGSetPrecond(solver, precond_id,
351 1 precond, ierr)
352
353
354 c Now setup and solve!
355 call HYPRE_ParCSRPCGSetup(solver, parcsr_A, par_b,
356 1 par_x, ierr)
357 call HYPRE_ParCSRPCGSolve(solver, parcsr_A, par_b,
358 1 par_x, ierr)
359
360
361 c Run info - needed logging turned on
362
363 call HYPRE_ParCSRPCGGetNumIterations(solver, num_iterations,
364 1 ierr)
365 call HYPRE_ParCSRPCGGetFinalRelative(solver, final_res_norm,
366 1 ierr)
367
368 if ( myid .eq. 0 ) then
369 print *
370 print *, "Iterations = ", num_iterations
371 print *, "Final Relative Residual Norm = ", final_res_norm
372 print *
373 endif
374
375 c Destroy precond and solver
376
377 call HYPRE_BoomerAMGDestroy(precond, ierr)
378 call HYPRE_ParCSRPCGDestroy(solver, ierr)
379
380 c PCG with ParaSails
381 elseif ( solver_id .eq. 8 ) then
382
383 c Create solver
384 call HYPRE_ParCSRPCGCreate(MPI_COMM_WORLD, solver, ierr)
385
386 c Set some parameters (See Reference Manual for more parameters)
387 call HYPRE_ParCSRPCGSetMaxIter(solver, 1000, ierr)
388 call HYPRE_ParCSRPCGSetTol(solver, 1.0d-7, ierr)
389 call HYPRE_ParCSRPCGSetTwoNorm(solver, 1, ierr)
390 call HYPRE_ParCSRPCGSetPrintLevel(solver, 2, ierr)
391 call HYPRE_ParCSRPCGSetLogging(solver, 1, ierr)
392
393 c Now set up the Parasails preconditioner and specify any parameters
394 call HYPRE_ParaSailsCreate(MPI_COMM_WORLD, precond,ierr)
395 call HYPRE_ParaSailsSetParams(precond, 0.1d0, 1, ierr)
396 call HYPRE_ParaSailsSetFilter(precond, 0.05d0, ierr)
397 call HYPRE_ParaSailsSetSym(precond, 1)
398 call HYPRE_ParaSailsSetLogging(precond, 3, ierr)
399
400 c set parsails as the pcg preconditioner
401 precond_id = 4
402 call HYPRE_ParCSRPCGSetPrecond(solver, precond_id,
403 1 precond, ierr)
404
405
406 c Now setup and solve!
407 call HYPRE_ParCSRPCGSetup(solver, parcsr_A, par_b,
408 1 par_x, ierr)
409 call HYPRE_ParCSRPCGSolve(solver, parcsr_A, par_b,
410 1 par_x, ierr)
411
412
413 c Run info - needed logging turned on
414
415 call HYPRE_ParCSRPCGGetNumIterations(solver, num_iterations,
416 1 ierr)
417 call HYPRE_ParCSRPCGGetFinalRelative(solver, final_res_norm,
418 1 ierr)
419
420 if ( myid .eq. 0 ) then
421 print *
422 print *, "Iterations = ", num_iterations
423 print *, "Final Relative Residual Norm = ", final_res_norm
424 print *
425 endif
426
427 c Destroy precond and solver
428
429 call HYPRE_ParaSailsDestroy(precond, ierr)
430 call HYPRE_ParCSRPCGDestroy(solver, ierr)
431
432 else
433 if ( myid .eq. 0 ) then
434 print *,'Invalid solver id specified'
435 stop
436 endif
437 endif
438
439
440
441 c Print the solution
442 if ( print_solution .ne. 0 ) then
443 call HYPRE_IJVectorPrint(x, "ij.out.x", ierr)
444 endif
445
446 c Clean up
447
448 call HYPRE_IJMatrixDestroy(A, ierr)
449 call HYPRE_IJVectorDestroy(b, ierr)
450 call HYPRE_IJVectorDestroy(x, ierr)
451
452
453 c Finalize MPI
454 call MPI_Finalize(ierr)
455
456 stop
457 end
syntax highlighted by Code2HTML, v. 0.9.1