download the original source code.
1 c
2 c Example 12
3 c
4 c Interface: Semi-Structured interface (SStruct)
5 c
6 c Compile with: make ex12f (may need to edit HYPRE_DIR in Makefile)
7 c
8 c Sample runs: mpirun -np 2 ex12f
9 c
10 c Description: The grid layout is the same as ex1, but with nodal
11 c unknowns. The solver is PCG preconditioned with either PFMG or
12 c BoomerAMG, set with 'precond_id' below.
13 c
14 c We recommend viewing the Struct examples before viewing this and
15 c the other SStruct examples. This is one of the simplest SStruct
16 c examples, used primarily to demonstrate how to set up
17 c non-cell-centered problems, and to demonstrate how easy it is to
18 c switch between structured solvers (PFMG) and solvers designed for
19 c more general settings (AMG).
20 c
21
22 program ex12f
23
24 implicit none
25
26 include 'mpif.h'
27 include 'HYPREf.h'
28
29 integer ierr
30 integer i, j, myid, num_procs
31
32 integer*8 grid
33 integer*8 graph
34 integer*8 stencil
35 integer*8 A
36 integer*8 b
37 integer*8 x
38
39 integer nparts
40 integer nvars
41 integer part
42 integer var
43
44 integer precond_id, object_type
45
46 integer ilower(2), iupper(2)
47 integer vartypes(1)
48 integer offsets(2,5)
49 integer ent
50 integer nentries, nvalues, stencil_indices(5)
51
52 double precision values(100), tol
53
54 c This comes from 'sstruct_mv/HYPRE_sstruct_mv.h'
55 integer HYPRE_SSTRUCT_VARIABLE_NODE
56 parameter( HYPRE_SSTRUCT_VARIABLE_NODE = 1 )
57
58 integer*8 sA
59 integer*8 sb
60 integer*8 sx
61 integer*8 parA
62 integer*8 parb
63 integer*8 parx
64 integer*8 solver
65 integer*8 precond
66
67 character*32 matfile
68
69 c We only have one part and one variable
70 nparts = 1
71 nvars = 1
72 part = 0
73 var = 0
74
75 c Initialize MPI
76 call MPI_Init(ierr)
77 call MPI_Comm_rank(MPI_COMM_WORLD, myid, ierr)
78 call MPI_Comm_size(MPI_COMM_WORLD, num_procs, ierr)
79
80 if (num_procs .ne. 2) then
81 if (myid .eq. 0) then
82 print *, "Must run with 2 processors!"
83 stop
84 endif
85 endif
86
87 c Set preconditioner id (PFMG = 1, BoomerAMG = 2)
88 precond_id = 1
89
90 if (precond_id .eq. 1) then
91 object_type = HYPRE_STRUCT
92 else if (precond_id .eq. 2) then
93 object_type = HYPRE_PARCSR
94 else
95 if (myid .eq. 0) then
96 print *, "Invalid solver!"
97 stop
98 endif
99 endif
100
101 c-----------------------------------------------------------------------
102 c 1. Set up the grid. Here we use only one part. Each processor
103 c describes the piece of the grid that it owns.
104 c-----------------------------------------------------------------------
105
106 c Create an empty 2D grid object
107 call HYPRE_SStructGridCreate(MPI_COMM_WORLD, 2, nparts, grid,
108 + ierr)
109
110 c Add boxes to the grid
111 if (myid .eq. 0) then
112 ilower(1) = -3
113 ilower(2) = 1
114 iupper(1) = -1
115 iupper(2) = 2
116 call HYPRE_SStructGridSetExtents(grid, part, ilower, iupper,
117 + ierr)
118 else if (myid .eq. 1) then
119 ilower(1) = 0
120 ilower(2) = 1
121 iupper(1) = 2
122 iupper(2) = 4
123 call HYPRE_SStructGridSetExtents(grid, part, ilower, iupper,
124 + ierr)
125 endif
126
127 c Set the variable type and number of variables on each part
128 vartypes(1) = HYPRE_SSTRUCT_VARIABLE_NODE
129 call HYPRE_SStructGridSetVariables(grid, part, nvars, vartypes,
130 + ierr)
131
132 c This is a collective call finalizing the grid assembly
133 call HYPRE_SStructGridAssemble(grid, ierr)
134
135 c-----------------------------------------------------------------------
136 c 2. Define the discretization stencil
137 c-----------------------------------------------------------------------
138
139 c Create an empty 2D, 5-pt stencil object
140 call HYPRE_SStructStencilCreate(2, 5, stencil, ierr)
141
142 c Define the geometry of the stencil. Each represents a relative
143 c offset (in the index space).
144 offsets(1,1) = 0
145 offsets(2,1) = 0
146 offsets(1,2) = -1
147 offsets(2,2) = 0
148 offsets(1,3) = 1
149 offsets(2,3) = 0
150 offsets(1,4) = 0
151 offsets(2,4) = -1
152 offsets(1,5) = 0
153 offsets(2,5) = 1
154
155 c Assign numerical values to the offsets so that we can easily refer
156 c to them - the last argument indicates the variable for which we
157 c are assigning this stencil
158 do ent = 1, 5
159 call HYPRE_SStructStencilSetEntry(stencil,
160 + ent-1, offsets(1,ent), var, ierr)
161 enddo
162
163 c-----------------------------------------------------------------------
164 c 3. Set up the Graph - this determines the non-zero structure of
165 c the matrix and allows non-stencil relationships between the parts
166 c-----------------------------------------------------------------------
167
168 c Create the graph object
169 call HYPRE_SStructGraphCreate(MPI_COMM_WORLD, grid, graph, ierr)
170
171 c See MatrixSetObjectType below
172 call HYPRE_SStructGraphSetObjectType(graph, object_type, ierr)
173
174 c Now we need to tell the graph which stencil to use for each
175 c variable on each part (we only have one variable and one part)
176 call HYPRE_SStructGraphSetStencil(graph, part, var, stencil, ierr)
177
178 c Here we could establish connections between parts if we had more
179 c than one part using the graph. For example, we could use
180 c HYPRE_GraphAddEntries() routine or HYPRE_GridSetNeighborPart()
181
182 c Assemble the graph
183 call HYPRE_SStructGraphAssemble(graph, ierr)
184
185 c-----------------------------------------------------------------------
186 c 4. Set up a SStruct Matrix
187 c-----------------------------------------------------------------------
188
189 c Create an empty matrix object
190 call HYPRE_SStructMatrixCreate(MPI_COMM_WORLD, graph, A, ierr)
191
192 c Set the object type (by default HYPRE_SSTRUCT). This determines
193 c the data structure used to store the matrix. For PFMG we use
194 c HYPRE_STRUCT, and for BoomerAMG we use HYPRE_PARCSR (set above).
195 call HYPRE_SStructMatrixSetObjectTyp(A, object_type, ierr)
196
197 c Get ready to set values
198 call HYPRE_SStructMatrixInitialize(A, ierr)
199
200 c Set the matrix coefficients. Each processor assigns coefficients
201 c for the boxes in the grid that it owns. Note that the
202 c coefficients associated with each stencil entry may vary from grid
203 c point to grid point if desired. Here, we first set the same
204 c stencil entries for each grid point. Then we make modifications
205 c to grid points near the boundary. Note that the ilower values are
206 c different from those used in ex1 because of the way nodal
207 c variables are referenced. Also note that some of the stencil
208 c values are set on both processor 0 and processor 1. See the User
209 c and Reference manuals for more details.
210
211 c Stencil entry labels correspond to the offsets defined above
212 do i = 1, 5
213 stencil_indices(i) = i-1
214 enddo
215 nentries = 5
216
217 if (myid .eq. 0) then
218 ilower(1) = -4
219 ilower(2) = 0
220 iupper(1) = -1
221 iupper(2) = 2
222 c 12 grid points, each with 5 stencil entries
223 nvalues = 60
224 else if (myid .eq. 1) then
225 ilower(1) = -1
226 ilower(2) = 0
227 iupper(1) = 2
228 iupper(2) = 4
229 c 12 grid points, each with 5 stencil entries
230 nvalues = 100
231 endif
232
233 do i = 1, nvalues, nentries
234 values(i) = 4.0
235 do j = 1, nentries-1
236 values(i+j) = -1.0
237 enddo
238 enddo
239
240 call HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
241 + var, nentries, stencil_indices, values, ierr)
242
243 c Set the coefficients reaching outside of the boundary to 0. Note
244 c that both ilower *and* iupper may be different from those in ex1.
245
246 do i = 1, 5
247 values(i) = 0.0
248 enddo
249
250 if (myid .eq. 0) then
251
252 c values below our box
253 ilower(1) = -4
254 ilower(2) = 0
255 iupper(1) = -1
256 iupper(2) = 0
257 stencil_indices(1) = 3
258 call HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
259 + var, 1, stencil_indices, values, ierr)
260 c values to the left of our box
261 ilower(1) = -4
262 ilower(2) = 0
263 iupper(1) = -4
264 iupper(2) = 2
265 stencil_indices(1) = 1
266 call HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
267 + var, 1, stencil_indices, values, ierr)
268 c values above our box
269 ilower(1) = -4
270 ilower(2) = 2
271 iupper(1) = -2
272 iupper(2) = 2
273 stencil_indices(1) = 4
274 call HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
275 + var, 1, stencil_indices, values, ierr)
276
277 else if (myid .eq. 1) then
278
279 c values below our box
280 ilower(1) = -1
281 ilower(2) = 0
282 iupper(1) = 2
283 iupper(2) = 0
284 stencil_indices(1) = 3
285 call HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
286 + var, 1, stencil_indices, values, ierr)
287 c values to the right of our box
288 ilower(1) = 2
289 ilower(2) = 0
290 iupper(1) = 2
291 iupper(2) = 4
292 stencil_indices(1) = 2
293 call HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
294 + var, 1, stencil_indices, values, ierr)
295 c values above our box
296 ilower(1) = -1
297 ilower(2) = 4
298 iupper(1) = 2
299 iupper(2) = 4
300 stencil_indices(1) = 4
301 call HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
302 + var, 1, stencil_indices, values, ierr)
303 c values to the left of our box
304 c (that do not border the other box on proc. 0)
305 ilower(1) = -1
306 ilower(2) = 3
307 iupper(1) = -1
308 iupper(2) = 4
309 stencil_indices(1) = 1
310 call HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
311 + var, 1, stencil_indices, values, ierr)
312
313 endif
314
315 c This is a collective call finalizing the matrix assembly
316 call HYPRE_SStructMatrixAssemble(A, ierr)
317
318 c matfile = 'ex12f.out'
319 c matfile(10:10) = char(0)
320 c call HYPRE_SStructMatrixPrint(matfile, A, 0, ierr)
321
322 c Create an empty vector object
323 call HYPRE_SStructVectorCreate(MPI_COMM_WORLD, grid, b, ierr)
324 call HYPRE_SStructVectorCreate(MPI_COMM_WORLD, grid, x, ierr)
325
326 c As with the matrix, set the appropriate object type for the vectors
327 call HYPRE_SStructVectorSetObjectTyp(b, object_type, ierr)
328 call HYPRE_SStructVectorSetObjectTyp(x, object_type, ierr)
329
330 c Indicate that the vector coefficients are ready to be set
331 call HYPRE_SStructVectorInitialize(b, ierr)
332 call HYPRE_SStructVectorInitialize(x, ierr)
333
334 c Set the vector coefficients. Again, note that the ilower values
335 c are different from those used in ex1, and some of the values are
336 c set on both processors.
337
338 if (myid .eq. 0) then
339
340 ilower(1) = -4
341 ilower(2) = 0
342 iupper(1) = -1
343 iupper(2) = 2
344
345 do i = 1, 12
346 values(i) = 1.0
347 enddo
348 call HYPRE_SStructVectorSetBoxValues(b, part, ilower, iupper,
349 + var, values, ierr)
350 do i = 1, 12
351 values(i) = 0.0
352 enddo
353 call HYPRE_SStructVectorSetBoxValues(x, part, ilower, iupper,
354 + var, values, ierr)
355
356 else if (myid .eq. 1) then
357
358 ilower(1) = 0
359 ilower(2) = 1
360 iupper(1) = 2
361 iupper(2) = 4
362
363 do i = 1, 20
364 values(i) = 1.0
365 enddo
366 call HYPRE_SStructVectorSetBoxValues(b, part, ilower, iupper,
367 + var, values, ierr)
368 do i = 1, 20
369 values(i) = 0.0
370 enddo
371 call HYPRE_SStructVectorSetBoxValues(x, part, ilower, iupper,
372 + var, values, ierr)
373
374 endif
375
376 c This is a collective call finalizing the vector assembly
377 call HYPRE_SStructVectorAssemble(b, ierr)
378 call HYPRE_SStructVectorAssemble(x, ierr)
379
380 c-----------------------------------------------------------------------
381 c 6. Set up and use a solver (See the Reference Manual for
382 c descriptions of all of the options.)
383 c-----------------------------------------------------------------------
384
385 tol = 1.0E-6
386
387 if (precond_id .eq. 1) then
388
389 c PFMG
390
391 c Because we are using a struct solver, we need to get the object
392 c of the matrix and vectors to pass in to the struct solvers
393 call HYPRE_SStructMatrixGetObject(A, sA, ierr)
394 call HYPRE_SStructVectorGetObject(b, sb, ierr)
395 call HYPRE_SStructVectorGetObject(x, sx, ierr)
396
397 c Create an empty PCG Struct solver
398 call HYPRE_StructPCGCreate(MPI_COMM_WORLD, solver, ierr)
399 c Set PCG parameters
400 call HYPRE_StructPCGSetTol(solver, tol, ierr)
401 call HYPRE_StructPCGSetPrintLevel(solver, 2, ierr)
402 call HYPRE_StructPCGSetMaxIter(solver, 50, ierr)
403
404 c Create the Struct PFMG solver for use as a preconditioner
405 call HYPRE_StructPFMGCreate(MPI_COMM_WORLD, precond, ierr)
406 c Set PFMG parameters
407 call HYPRE_StructPFMGSetMaxIter(precond, 1, ierr)
408 call HYPRE_StructPFMGSetTol(precond, 0.0, ierr)
409 call HYPRE_StructPFMGSetZeroGuess(precond, ierr)
410 call HYPRE_StructPFMGSetNumPreRelax(precond, 2, ierr)
411 call HYPRE_StructPFMGSetNumPostRelax(precond, 2, ierr)
412 c Non-Galerkin coarse grid (more efficient for this problem)
413 call HYPRE_StructPFMGSetRAPType(precond, 1, ierr)
414 c R/B Gauss-Seidel
415 call HYPRE_StructPFMGSetRelaxType(precond, 2, ierr)
416 c Skip relaxation on some levels (more efficient for this problem)
417 call HYPRE_StructPFMGSetSkipRelax(precond, 1, ierr)
418 c Set preconditioner (PFMG = 1) and solve
419 call HYPRE_StructPCGSetPrecond(solver, 1, precond, ierr)
420 call HYPRE_StructPCGSetup(solver, sA, sb, sx, ierr)
421 call HYPRE_StructPCGSolve(solver, sA, sb, sx, ierr)
422
423 c Free memory
424 call HYPRE_StructPCGDestroy(solver, ierr)
425 call HYPRE_StructPFMGDestroy(precond, ierr)
426
427 else if (precond_id .eq. 2) then
428
429 c BoomerAMG
430
431 c Because we are using a struct solver, we need to get the object
432 c of the matrix and vectors to pass in to the struct solvers
433 call HYPRE_SStructMatrixGetObject(A, parA, ierr)
434 call HYPRE_SStructVectorGetObject(b, parb, ierr)
435 call HYPRE_SStructVectorGetObject(x, parx, ierr)
436
437 c Create an empty PCG Struct solver
438 call HYPRE_ParCSRPCGCreate(MPI_COMM_WORLD, solver, ierr)
439 c Set PCG parameters
440 call HYPRE_ParCSRPCGSetTol(solver, tol, ierr)
441 call HYPRE_ParCSRPCGSetPrintLevel(solver, 2, ierr)
442 call HYPRE_ParCSRPCGSetMaxIter(solver, 50, ierr)
443
444 c Create the BoomerAMG solver for use as a preconditioner
445 call HYPRE_BoomerAMGCreate(precond, ierr)
446 c Set BoomerAMG parameters
447 call HYPRE_BoomerAMGSetMaxIter(precond, 1, ierr)
448 call HYPRE_BoomerAMGSetTol(precond, 0.0, ierr)
449 c Print amg solution info
450 call HYPRE_BoomerAMGSetPrintLevel(precond, 1, ierr)
451 call HYPRE_BoomerAMGSetCoarsenType(precond, 6, ierr)
452 call HYPRE_BoomerAMGSetOldDefault(precond, ierr)
453 c Sym G.S./Jacobi hybrid
454 call HYPRE_BoomerAMGSetRelaxType(precond, 6, ierr)
455 call HYPRE_BoomerAMGSetNumSweeps(precond, 1, ierr)
456 c Set preconditioner (BoomerAMG = 2) and solve
457 call HYPRE_ParCSRPCGSetPrecond(solver, 2, precond, ierr)
458 call HYPRE_ParCSRPCGSetup(solver, parA, parb, parx, ierr)
459 call HYPRE_ParCSRPCGSolve(solver, parA, parb, parx, ierr)
460
461 c Free memory
462 call HYPRE_ParCSRPCGDestroy(solver, ierr)
463 call HYPRE_BoomerAMGDestroy(precond, ierr)
464
465 endif
466
467 c Free memory
468 call HYPRE_SStructGridDestroy(grid, ierr)
469 call HYPRE_SStructStencilDestroy(stencil, ierr)
470 call HYPRE_SStructGraphDestroy(graph, ierr)
471 call HYPRE_SStructMatrixDestroy(A, ierr)
472 call HYPRE_SStructVectorDestroy(b, ierr)
473 call HYPRE_SStructVectorDestroy(x, ierr)
474
475 c Finalize MPI
476 call MPI_Finalize(ierr)
477
478 stop
479 end
syntax highlighted by Code2HTML, v. 0.9.1