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