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 ex12 (may need to edit HYPRE_DIR in Makefile)
7 c
8 c Sample runs: mpirun -np 2 ex12
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 Now we need to tell the graph which stencil to use for each
172 c variable on each part (we only have one variable and one part)
173 call HYPRE_SStructGraphSetStencil(graph, part, var, stencil, ierr)
174
175 c Here we could establish connections between parts if we had more
176 c than one part using the graph. For example, we could use
177 c HYPRE_GraphAddEntries() routine or HYPRE_GridSetNeighborPart()
178
179 c Assemble the graph
180 call HYPRE_SStructGraphAssemble(graph, ierr)
181
182 c-----------------------------------------------------------------------
183 c 4. Set up a SStruct Matrix
184 c-----------------------------------------------------------------------
185
186 c Create an empty matrix object
187 call HYPRE_SStructMatrixCreate(MPI_COMM_WORLD, graph, A, ierr)
188
189 c Set the object type (by default HYPRE_SSTRUCT). This determines
190 c the data structure used to store the matrix. For PFMG we use
191 c HYPRE_STRUCT, and for BoomerAMG we use HYPRE_PARCSR (set above).
192 call HYPRE_SStructMatrixSetObjectTyp(A, object_type, ierr)
193
194 c Get ready to set values
195 call HYPRE_SStructMatrixInitialize(A, ierr)
196
197 c Set the matrix coefficients. Each processor assigns coefficients
198 c for the boxes in the grid that it owns. Note that the
199 c coefficients associated with each stencil entry may vary from grid
200 c point to grid point if desired. Here, we first set the same
201 c stencil entries for each grid point. Then we make modifications
202 c to grid points near the boundary. Note that the ilower values are
203 c different from those used in ex1 because of the way nodal
204 c variables are referenced. Also note that some of the stencil
205 c values are set on both processor 0 and processor 1. See the User
206 c and Reference manuals for more details.
207
208 c Stencil entry labels correspond to the offsets defined above
209 do i = 1, 5
210 stencil_indices(i) = i-1
211 enddo
212 nentries = 5
213
214 if (myid .eq. 0) then
215 ilower(1) = -4
216 ilower(2) = 0
217 iupper(1) = -1
218 iupper(2) = 2
219 nvalues = 60 /* 12 grid points, each with 5 stencil entries */
220 else if (myid .eq. 1) then
221 ilower(1) = -1
222 ilower(2) = 0
223 iupper(1) = 2
224 iupper(2) = 4
225 nvalues = 100 /* 12 grid points, each with 5 stencil entries */
226 endif
227
228 do i = 1, nvalues, nentries
229 values(i) = 4.0
230 do j = 1, nentries-1
231 values(i+j) = -1.0
232 enddo
233 enddo
234
235 call HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
236 + var, nentries, stencil_indices, values, ierr)
237
238 c Set the coefficients reaching outside of the boundary to 0. Note
239 c that both ilower *and* iupper may be different from those in ex1.
240
241 do i = 1, 5
242 values(i) = 0.0
243 enddo
244
245 if (myid .eq. 0) then
246
247 c values below our box
248 ilower(1) = -4
249 ilower(2) = 0
250 iupper(1) = -1
251 iupper(2) = 0
252 stencil_indices(1) = 3
253 call HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
254 + var, 1, stencil_indices, values, ierr)
255 c values to the left of our box
256 ilower(1) = -4
257 ilower(2) = 0
258 iupper(1) = -4
259 iupper(2) = 2
260 stencil_indices(1) = 1
261 call HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
262 + var, 1, stencil_indices, values, ierr)
263 c values above our box
264 ilower(1) = -4
265 ilower(2) = 2
266 iupper(1) = -2
267 iupper(2) = 2
268 stencil_indices(1) = 4
269 call HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
270 + var, 1, stencil_indices, values, ierr)
271
272 else if (myid .eq. 1) then
273
274 c values below our box
275 ilower(1) = -1
276 ilower(2) = 0
277 iupper(1) = 2
278 iupper(2) = 0
279 stencil_indices(1) = 3
280 call HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
281 + var, 1, stencil_indices, values, ierr)
282 c values to the right of our box
283 ilower(1) = 2
284 ilower(2) = 0
285 iupper(1) = 2
286 iupper(2) = 4
287 stencil_indices(1) = 2
288 call HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
289 + var, 1, stencil_indices, values, ierr)
290 c values above our box
291 ilower(1) = -1
292 ilower(2) = 4
293 iupper(1) = 2
294 iupper(2) = 4
295 stencil_indices(1) = 4
296 call HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
297 + var, 1, stencil_indices, values, ierr)
298 c values to the left of our box
299 c (that do not border the other box on proc. 0)
300 ilower(1) = -1
301 ilower(2) = 3
302 iupper(1) = -1
303 iupper(2) = 4
304 stencil_indices(1) = 1
305 call HYPRE_SStructMatrixSetBoxValues(A, part, ilower, iupper,
306 + var, 1, stencil_indices, values, ierr)
307
308 endif
309
310 c This is a collective call finalizing the matrix assembly
311 call HYPRE_SStructMatrixAssemble(A, ierr)
312
313 c matfile = 'ex12f.out'
314 c matfile(10:10) = char(0)
315 c call HYPRE_SStructMatrixPrint(matfile, A, 0, ierr)
316
317 c Create an empty vector object
318 call HYPRE_SStructVectorCreate(MPI_COMM_WORLD, grid, b, ierr)
319 call HYPRE_SStructVectorCreate(MPI_COMM_WORLD, grid, x, ierr)
320
321 c As with the matrix, set the appropriate object type for the vectors
322 call HYPRE_SStructVectorSetObjectTyp(b, object_type, ierr)
323 call HYPRE_SStructVectorSetObjectTyp(x, object_type, ierr)
324
325 c Indicate that the vector coefficients are ready to be set
326 call HYPRE_SStructVectorInitialize(b, ierr)
327 call HYPRE_SStructVectorInitialize(x, ierr)
328
329 c Set the vector coefficients. Again, note that the ilower values
330 c are different from those used in ex1, and some of the values are
331 c set on both processors.
332
333 if (myid .eq. 0) then
334
335 ilower(1) = -4
336 ilower(2) = 0
337 iupper(1) = -1
338 iupper(2) = 2
339
340 do i = 1, 12
341 values(i) = 1.0
342 enddo
343 call HYPRE_SStructVectorSetBoxValues(b, part, ilower, iupper,
344 + var, values, ierr)
345 do i = 1, 12
346 values(i) = 0.0
347 enddo
348 call HYPRE_SStructVectorSetBoxValues(x, part, ilower, iupper,
349 + var, values, ierr)
350
351 else if (myid .eq. 1) then
352
353 ilower(1) = 0
354 ilower(2) = 1
355 iupper(1) = 2
356 iupper(2) = 4
357
358 do i = 1, 20
359 values(i) = 1.0
360 enddo
361 call HYPRE_SStructVectorSetBoxValues(b, part, ilower, iupper,
362 + var, values, ierr)
363 do i = 1, 20
364 values(i) = 0.0
365 enddo
366 call HYPRE_SStructVectorSetBoxValues(x, part, ilower, iupper,
367 + var, values, ierr)
368
369 endif
370
371 c This is a collective call finalizing the vector assembly
372 call HYPRE_SStructVectorAssemble(b, ierr)
373 call HYPRE_SStructVectorAssemble(x, ierr)
374
375 c-----------------------------------------------------------------------
376 c 6. Set up and use a solver (See the Reference Manual for
377 c descriptions of all of the options.)
378 c-----------------------------------------------------------------------
379
380 tol = 1.0E-6
381
382 if (precond_id .eq. 1) then
383
384 c PFMG
385
386 c Because we are using a struct solver, we need to get the object
387 c of the matrix and vectors to pass in to the struct solvers
388 call HYPRE_SStructMatrixGetObject(A, sA, ierr)
389 call HYPRE_SStructVectorGetObject(b, sb, ierr)
390 call HYPRE_SStructVectorGetObject(x, sx, ierr)
391
392 c Create an empty PCG Struct solver
393 call HYPRE_StructPCGCreate(MPI_COMM_WORLD, solver, ierr)
394 c Set PCG parameters
395 call HYPRE_StructPCGSetTol(solver, tol, ierr)
396 call HYPRE_StructPCGSetPrintLevel(solver, 2, ierr)
397 call HYPRE_StructPCGSetMaxIter(solver, 50, ierr)
398
399 c Create the Struct PFMG solver for use as a preconditioner
400 call HYPRE_StructPFMGCreate(MPI_COMM_WORLD, precond, ierr)
401 c Set PFMG parameters
402 call HYPRE_StructPFMGSetMaxIter(precond, 1, ierr)
403 call HYPRE_StructPFMGSetTol(precond, 0.0, ierr)
404 call HYPRE_StructPFMGSetZeroGuess(precond, ierr)
405 call HYPRE_StructPFMGSetNumPreRelax(precond, 2, ierr)
406 call HYPRE_StructPFMGSetNumPostRelax(precond, 2, ierr)
407 c Non-Galerkin coarse grid (more efficient for this problem)
408 call HYPRE_StructPFMGSetRAPType(precond, 1, ierr)
409 c R/B Gauss-Seidel
410 call HYPRE_StructPFMGSetRelaxType(precond, 2, ierr)
411 c Skip relaxation on some levels (more efficient for this problem)
412 call HYPRE_StructPFMGSetSkipRelax(precond, 1, ierr)
413 c Set preconditioner (PFMG = 1) and solve
414 call HYPRE_StructPCGSetPrecond(solver, 1, precond, ierr)
415 call HYPRE_StructPCGSetup(solver, sA, sb, sx, ierr)
416 call HYPRE_StructPCGSolve(solver, sA, sb, sx, ierr)
417
418 c Free memory
419 call HYPRE_StructPCGDestroy(solver, ierr)
420 call HYPRE_StructPFMGDestroy(precond, ierr)
421
422 else if (precond_id .eq. 2) then
423
424 c BoomerAMG
425
426 c Because we are using a struct solver, we need to get the object
427 c of the matrix and vectors to pass in to the struct solvers
428 call HYPRE_SStructMatrixGetObject(A, parA, ierr)
429 call HYPRE_SStructVectorGetObject(b, parb, ierr)
430 call HYPRE_SStructVectorGetObject(x, parx, ierr)
431
432 c Create an empty PCG Struct solver
433 call HYPRE_ParCSRPCGCreate(MPI_COMM_WORLD, solver, ierr)
434 c Set PCG parameters
435 call HYPRE_ParCSRPCGSetTol(solver, tol, ierr)
436 call HYPRE_ParCSRPCGSetPrintLevel(solver, 2, ierr)
437 call HYPRE_ParCSRPCGSetMaxIter(solver, 50, ierr)
438
439 c Create the BoomerAMG solver for use as a preconditioner
440 call HYPRE_BoomerAMGCreate(precond, ierr)
441 c Set BoomerAMG parameters
442 call HYPRE_BoomerAMGSetMaxIter(precond, 1, ierr)
443 call HYPRE_BoomerAMGSetTol(precond, 0.0, ierr)
444 c Print amg solution info
445 call HYPRE_BoomerAMGSetPrintLevel(precond, 1, ierr)
446 call HYPRE_BoomerAMGSetCoarsenType(precond, 6, ierr)
447 c Sym G.S./Jacobi hybrid
448 call HYPRE_BoomerAMGSetRelaxType(precond, 6, ierr)
449 call HYPRE_BoomerAMGSetNumSweeps(precond, 1, ierr)
450 c Set preconditioner (BoomerAMG = 2) and solve
451 call HYPRE_ParCSRPCGSetPrecond(solver, 2, precond, ierr)
452 call HYPRE_ParCSRPCGSetup(solver, parA, parb, parx, ierr)
453 call HYPRE_ParCSRPCGSolve(solver, parA, parb, parx, ierr)
454
455 c Free memory
456 call HYPRE_ParCSRPCGDestroy(solver, ierr)
457 call HYPRE_BoomerAMGDestroy(precond, ierr)
458
459 endif
460
461 c Free memory
462 call HYPRE_SStructGridDestroy(grid, ierr)
463 call HYPRE_SStructStencilDestroy(stencil, ierr)
464 call HYPRE_SStructGraphDestroy(graph, ierr)
465 call HYPRE_SStructMatrixDestroy(A, ierr)
466 call HYPRE_SStructVectorDestroy(b, ierr)
467 call HYPRE_SStructVectorDestroy(x, ierr)
468
469 c Finalize MPI
470 call MPI_Finalize(ierr)
471
472 stop
473 end
syntax highlighted by Code2HTML, v. 0.9.1