-
Notifications
You must be signed in to change notification settings - Fork 87
/
valout.f90
499 lines (423 loc) · 20.2 KB
/
valout.f90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
!! Geoclaw specific output - adds eta to q array before writing out
!!
!! Write the results to the file fort.q<iframe>
!! Use format required by matlab script plotclaw2.m or Python tools
!!
!! set outaux = .true. to also output the aux arrays to fort.a<iframe>
subroutine valout(level_begin, level_end, time, num_eqn, num_aux)
use amr_module, only: alloc, t0, output_aux_onlyonce, output_aux_components
use amr_module, only: frame => matlabu, num_ghost => nghost, lstart
use amr_module, only: hxposs, hyposs, output_format, store1, storeaux
use amr_module, only: node, rnode, ndilo, ndihi, ndjlo, ndjhi
use amr_module, only: cornxlo, cornylo, levelptr, mxnest
use amr_module, only: timeValout, timeValoutCPU, tvoll, tvollCPU, rvoll
use amr_module, only: timeTick, tick_clock_start, t0, timeTickCPU
use storm_module, only: storm_specification_type, output_storm_location
use storm_module, only: output_storm_location
use storm_module, only: landfall, display_landfall_time
#ifdef HDF5
use hdf5
#endif
implicit none
! Input
integer, intent(in) :: level_begin, level_end, num_eqn, num_aux
real(kind=8), intent(in) :: time
! Locals
logical :: timing_file_exists
integer, parameter :: out_unit = 50
integer, parameter :: binary_unit = 51
integer :: i, j, m, level, output_aux_num, num_stop, digit
integer :: grid_ptr, num_cells(2), num_grids, q_loc, aux_loc
real(kind=8) :: lower_corner(2), delta(2)
logical :: out_aux
character(len=11) :: file_name(5)
real(kind=8) :: h, hu, hv, eta
real(kind=8), allocatable :: qeta(:)
real(kind=4), allocatable :: qeta4(:), aux4(:)
integer :: lenaux4
#ifdef HDF5
! HDF writing
integer :: hdf_error
integer(hid_t) :: hdf_file, data_space, data_set
integer(hsize_t) :: dims(2)
#endif
! Timing
integer(kind=8) :: clock_start, clock_finish, clock_rate
integer(kind=8) :: tick_clock_finish, tick_clock_rate, timeTick_int
real(kind=8) :: cpu_start, cpu_finish, t_CPU_overall, timeTick_overall
character(len=128) :: console_format
character(len=512) :: timing_line, timing_substr
character(len=*), parameter :: timing_file_name = "timing.csv"
character(len=8) :: file_format
character(len=*), parameter :: header_format = &
"(i6,' grid_number',/," // &
"i6,' AMR_level',/," // &
"i6,' mx',/," // &
"i6,' my',/" // &
"e26.16,' xlow', /, " // &
"e26.16,' ylow', /," // &
"e26.16,' dx', /," // &
"e26.16,' dy',/)"
character(len=*), parameter :: t_file_format = "(e18.8,' time', /," // &
"i6,' meqn'/," // &
"i6,' ngrids'/," // &
"i6,' naux'/," // &
"i6,' ndim'/," // &
"i6,' nghost'/," // &
"a10,' format'/,/)"
! Output timing
call system_clock(clock_start,clock_rate)
call cpu_time(cpu_start)
! Count how many aux components requested
output_aux_num = 0
do i=1, num_aux
output_aux_num = output_aux_num + output_aux_components(i)
end do
! Note: Currently outputs all aux components if any are requested
out_aux = ((output_aux_num > 0) .and. &
((.not. output_aux_onlyonce) .or. (abs(time - t0) < 1d-90)))
! Output storm track if needed
if (storm_specification_type /= 0) then
call output_storm_location(time)
end if
! Construct file names
file_name(1) = 'fort.qxxxx'
file_name(2) = 'fort.txxxx'
file_name(3) = 'fort.axxxx'
file_name(4) = 'fort.bxxxx'
num_stop = frame
do i = 10, 7, -1
digit = mod(num_stop, 10)
do j = 1, 4
file_name(j)(i:i) = char(ichar('0') + digit)
end do
num_stop = num_stop / 10
end do
! Slightly modified for HDF file output
file_name(5) = 'clawxxxx.h5'
num_stop = frame
do i = 8, 5, -1
digit = mod(num_stop, 10)
file_name(5)(i:i) = char(ichar('0') + digit)
num_stop = num_stop / 10
end do
! ==========================================================================
! Write out fort.q file (and fort.bXXXX and clawxxxx.h5 files if necessary)
! Here we let fort.q be out_unit and the the other two be out_unit + 1
open(unit=out_unit, file=file_name(1), status='unknown', form='formatted')
if (output_format == 2 .or. output_format == 3) then
open(unit=binary_unit, file=file_name(4), status="unknown", &
access='stream')
#ifdef HDF5
else if (output_format == 4) then
! Note that we will use this file for both q and aux data
call h5create_f(file_name(5), H5F_ACC_TRUNC_F, hdf_file, hdf_error)
! Create group for q
call h5gcreate_f(hdf_file, "/q", q_group, hdf_error)
#endif
end if
num_grids = 0
! Loop over levels
do level = level_begin, level_end
grid_ptr = lstart(level)
delta = [hxposs(level), hyposs(level)]
! Loop over grids on each level
do while (grid_ptr /= 0)
! Extract grid data
num_grids = num_grids + 1
num_cells(1) = node(ndihi, grid_ptr) - node(ndilo, grid_ptr) + 1
num_cells(2) = node(ndjhi, grid_ptr) - node(ndjlo, grid_ptr) + 1
q_loc = node(store1, grid_ptr)
aux_loc = node(storeaux, grid_ptr)
lower_corner = [rnode(cornxlo, grid_ptr), rnode(cornylo, grid_ptr)]
! Write out header data
write(out_unit, header_format) grid_ptr, level, &
num_cells(1), &
num_cells(2), &
lower_corner(1), &
lower_corner(2), &
delta(1), delta(2)
! Output grids
! Round off if nearly zero
! (Do this for all output_format's)
forall (m = 1:num_eqn, &
i=num_ghost + 1:num_cells(1) + num_ghost, &
j=num_ghost + 1:num_cells(2) + num_ghost, &
abs(alloc(iadd(m, i, j))) < 1d-90)
alloc(iadd(m, i, j)) = 0.d0
end forall
if (output_format == 1) then
! ascii output
do j = num_ghost + 1, num_cells(2) + num_ghost
do i = num_ghost + 1, num_cells(1) + num_ghost
! Extract depth and momenta
h = alloc(iadd(1, i, j))
hu = alloc(iadd(2, i, j))
hv = alloc(iadd(3, i, j))
! Calculate sufaces
eta = h + alloc(iaddaux(1,i,j))
if (abs(eta) < 1d-99) then
eta = 0.d0
end if
write(out_unit, "(50e26.16)") h, hu, hv, eta
end do
write(out_unit, *) ' '
end do
else if (output_format==2 .or. output_format==3) then
! binary32 or binary64
! Note: We are writing out ghost cell data also,
! so need to update this
call bound(time,num_eqn,num_ghost,alloc(q_loc), &
num_cells(1) + 2*num_ghost, &
num_cells(2) + 2*num_ghost, &
grid_ptr,alloc(aux_loc),num_aux)
! Need to add eta to the output data
allocate(qeta((num_eqn + 1) &
* (num_cells(1) + 2 * num_ghost) &
* (num_cells(2) + 2 * num_ghost)))
do j = 1, num_cells(2) + 2 * num_ghost
do i = 1, num_cells(1) + 2 * num_ghost
do m = 1, num_eqn
qeta(iaddqeta(m, i, j)) = alloc(iadd(m, i, j))
end do
eta = alloc(iadd(1, i, j)) + alloc(iaddaux(1, i ,j))
qeta(iaddqeta(num_eqn + 1, i, j)) = eta
end do
end do
if (output_format==2) then
! binary32 (shorten to 4-byte)
allocate(qeta4((num_eqn + 1) &
* (num_cells(1) + 2 * num_ghost) &
* (num_cells(2) + 2 * num_ghost)))
qeta4 = real(qeta, kind=4)
write(binary_unit) qeta4
deallocate(qeta4)
else if (output_format==3) then
! binary64 (full 8-byte)
write(binary_unit) qeta
endif
deallocate(qeta)
else if (output_format == 4) then
! HDF5 output
#ifdef HDF5
! Create data space - handles dimensions of the corresponding
! data set - annoyingling need to stick grid size into other
! data type
dims = (/ num_eqn, num_cells(1) + 2 * num_ghost, &
num_cells(2) + 2 * num_ghost /)
call h5screate_simple_f(2, dims, data_space, hdf_error)
! Create new dataset for this grid
call h5dcreate_f(hdf_file, data_set, H5T_IEEE_F64LE, &
data_space, data_set, hdf_error)
! Write q into file
i = (iadd(num_eqn, num_cells(1) + 2 * num_ghost, &
num_cells(2) + 2 * num_ghost))
call h5dwrite_f(data_set, H5T_NATIVE_DOUBLE, &
alloc(iadd(1, 1, 1):i), hdf_error)
call h5dclose_f(data_set, hdf_error)
call h5sclose_f(data_space, hdf_error)
#endif
else
print *, "Unsupported output format", output_format,"."
stop
endif
grid_ptr = node(levelptr, grid_ptr)
end do
end do
close(out_unit)
if (output_format==2 .or. output_format==3) then
close(binary_unit)
#ifdef HDF5
else if (output_format == 4) then
call h5gclose_f(q_group, hdf_error)
#endif
end if
! ==========================================================================
! Write out fort.a file
if (out_aux) then
if (output_format == 1) then
open(unit=out_unit, file=file_name(3), status='unknown', &
form='formatted')
else if (output_format==2 .or. output_format==3) then
open(unit=out_unit, file=file_name(3), status='unknown', &
access='stream')
#ifdef HDF5
else if (output_format == 4) then
! Create group for aux
call h5gcreate_f(hdf_file, "/aux", aux_group, hdf_error)
#endif
end if
do level = level_begin, level_end
grid_ptr = lstart(level)
delta = [hxposs(level), hyposs(level)]
! Loop over grids on each level
do while (grid_ptr /= 0)
! Extract grid data
num_cells(1) = node(ndihi, grid_ptr) - node(ndilo, grid_ptr) + 1
num_cells(2) = node(ndjhi, grid_ptr) - node(ndjlo, grid_ptr) + 1
aux_loc = node(storeaux, grid_ptr)
lower_corner = [rnode(cornxlo, grid_ptr), &
rnode(cornylo, grid_ptr)]
! Output grids
select case(output_format)
! ASCII output
case(1)
! We only output header info for aux data if writing
! ASCII data
write(out_unit, header_format) grid_ptr, level, &
num_cells(1), &
num_cells(2), &
lower_corner(1), &
lower_corner(2), &
delta(1), delta(2)
! Round off if nearly zero
forall (m = 1:num_aux, &
i=num_ghost + 1:num_cells(1) + num_ghost, &
j=num_ghost + 1:num_cells(2) + num_ghost, &
abs(alloc(iaddaux(m, i, j))) < 1d-90)
alloc(iaddaux(m, i, j)) = 0.d0
end forall
do j = num_ghost + 1, num_cells(2) + num_ghost
do i = num_ghost + 1, num_cells(1) + num_ghost
write(out_unit, "(50e26.16)") &
(alloc(iaddaux(m, i, j)), m=1, num_aux)
end do
write(out_unit, *) ' '
end do
case(2)
! binary32
! Note: We are writing out ghost cell data also
i = (iaddaux(num_aux, num_cells(1) + 2 * num_ghost, &
num_cells(2) + 2 * num_ghost))
lenaux4 = i - iaddaux(1, 1, 1) + 1
allocate(aux4(lenaux4))
aux4 = real(alloc(iaddaux(1, 1, 1):i), kind=4)
write(out_unit) aux4
deallocate(aux4)
case(3)
! binary64
! Note: We are writing out ghost cell data also
i = (iaddaux(num_aux, num_cells(1) + 2 * num_ghost, &
num_cells(2) + 2 * num_ghost))
write(out_unit) alloc(iaddaux(1, 1, 1):i)
! HDF5 output
case(4)
#ifdef HDF5
! Create data space - handles dimensions of the corresponding
! data set - annoyingling need to stick grid size into other
! data type
dims = (/ num_aux, num_cells(1) + 2 * num_ghost, &
num_cells(2) + 2 * num_ghost /)
call h5screate_simple_f(2, dims, data_space, hdf_error)
! Create new dataset for this grid
call h5dcreate_f(hdf_file, data_set, H5T_IEEE_F64LE, &
data_space, data_set, hdf_error)
call h5dcreate_f(hdf_file, data_set, H5T_IEEE_F64LE, &
data_space, data_set, hdf_error)
! Write q into file
i = (iadd_aux(num_aux, num_cells(1) + 2 * num_ghost, &
num_cells(2) + 2 * num_ghost))
call h5dwrite_f(data_set, H5T_NATIVE_DOUBLE, &
alloc(iadd_aux(1, 1, 1):i), hdf_error)
call h5dclose_f(data_set, hdf_error)
call h5sclose_f(data_space, hdf_error)
#endif
case default
print *, "Unsupported output format", output_format,"."
stop
end select
grid_ptr = node(levelptr, grid_ptr)
end do
end do
end if
close(out_unit)
#ifdef HDF5
if (out_aux) then
call h5gclose_f(aux_group, hdf_error)
end if
call h5fclose_f(hdf_file, hdf_error)
#endif
! ==========================================================================
! Write fort.t file
open(unit=out_unit, file=file_name(2), status='unknown', form='formatted')
! Note: We need to print out num_ghost too in order to strip ghost cells
! from q array when reading in pyclaw.io.binary
if (output_format == 1) then
file_format = 'ascii'
else if (output_format == 2) then
file_format = 'binary32'
else if (output_format == 3) then
file_format = 'binary64'
else if (output_format == 4) then
file_format = 'hdf'
endif
write(out_unit, t_file_format) time, num_eqn + 1, num_grids, num_aux, &
2, num_ghost, file_format
close(out_unit)
! ==========================================================================
! Write out timing stats
open(unit=out_unit, file=timing_file_name, form='formatted', &
status='unknown', action='write', position='append')
!status='old', action='write', position='append')
timing_line = "(e16.6, ', ', e16.6, ', ', e16.6,"
do level=1, mxnest
timing_substr = "', ', e16.6, ', ', e16.6, ', ', e16.6"
timing_line = trim(timing_line) // timing_substr
end do
timing_line = trim(timing_line) // ")"
if (abs(time - t0) < 1d-15) then
t_CPU_overall = 0.d0
timeTick_overall = 0.d0
else
call cpu_time(t_CPU_overall)
! if this is a restart, need to adjust add in time from previous run:
t_CPU_overall = t_CPU_overall + timeTickCPU
call system_clock(tick_clock_finish,tick_clock_rate)
timeTick_int = timeTick + tick_clock_finish - tick_clock_start
timeTick_overall = real(timeTick_int, kind=8)/real(clock_rate,kind=8)
endif
write(out_unit, timing_line) time, timeTick_overall, t_CPU_overall, &
(real(tvoll(i), kind=8) / real(clock_rate, kind=8), &
tvollCPU(i), rvoll(i), i=1,mxnest)
close(out_unit)
! ==========================================================================
! Print output info
if (display_landfall_time) then
! Convert time to days relative to landfall
console_format = "('GEOCLAW: Frame ',i4,' output files done at " // &
"time t = ', f5.2,/)"
print console_format, frame, time / (24.d0 * 60d0**2)
else
console_format = "('GEOCLAW: Frame ',i4,' output files done at " // &
"time t = ', d13.6,/)"
print console_format, frame, time
end if
! Increment frame counter
frame = frame + 1
! Ouptut timing
call system_clock(clock_finish,clock_rate)
call cpu_time(cpu_finish)
timeValout = timeValout + clock_finish - clock_start
timeValoutCPU = timeValoutCPU + cpu_finish - cpu_start
contains
! Index into q array
pure integer function iadd(m, i, j)
implicit none
integer, intent(in) :: m, i, j
iadd = q_loc + m - 1 + num_eqn * ((j - 1) * (num_cells(1) + 2 * num_ghost) + i - 1)
end function iadd
! Index into aux array
pure integer function iaddaux(m, i, j)
implicit none
integer, intent(in) :: m, i, j
iaddaux = aux_loc + m - 1 + num_aux * (i - 1) + num_aux * (num_cells(1) + 2 * num_ghost) * (j - 1)
end function iaddaux
! Index into qeta for binary output
! Note that this implicitly assumes that we are outputting only h, hu, hv
! and will not output more (change num_eqn parameter above)
pure integer function iaddqeta(m, i, j)
implicit none
integer, intent(in) :: m, i, j
iaddqeta = 1 + m - 1 + (num_eqn + 1) * ((j - 1) * (num_cells(1) + 2 * num_ghost) + i - 1)
end function iaddqeta
end subroutine valout