44 logical,
private :: module_setup = .false.
46 integer,
parameter :: OUTGAUGEUNIT = 89
49 integer,
parameter :: MAX_BUFFER = 1000
56 character(len=14) :: file_name
59 real(kind=8) :: x, y, t_start, t_end
62 real(kind=8) :: last_time
65 integer :: file_format
66 real(kind=8) :: min_time_increment
67 character(len=10) :: display_format
68 logical,
allocatable :: q_out_vars(:)
69 logical,
allocatable :: aux_out_vars(:)
70 integer :: num_out_vars
73 real(kind=8),
allocatable :: data(:, :)
74 integer :: level(MAX_BUFFER)
77 integer :: buffer_index
84 integer,
allocatable,
dimension(:) :: mbestsrc, mbestorder, &
85 igauge, mbestg1, mbestg2
92 use utility_module
, only: get_value_count
97 logical,
intent(in) :: restart
98 integer :: num_eqn, num_aux
99 character(len=*),
intent(in),
optional :: fname
102 integer :: i, n, index
103 integer :: num, pos, digit
104 integer,
parameter :: unit = 7
105 character(len=128) :: header_1
106 character(len=40) :: q_column, aux_column
108 if (.not. module_setup)
then
111 if (present(fname))
then
117 read(unit, *) num_gauges
118 allocate(gauges(num_gauges))
121 allocate(mbestsrc(num_gauges), mbestorder(num_gauges))
122 allocate(mbestg1(maxgr), mbestg2(maxgr))
127 read(unit, *) gauges(i)%gauge_num, gauges(i)%x, gauges(i)%y, &
128 gauges(i)%t_start, gauges(i)%t_end
129 gauges(i)%buffer_index = 1
130 gauges(i)%last_time = gauges(i)%t_start
136 read(unit, *) (gauges(i)%file_format, i=1, num_gauges)
139 read(unit, *) (gauges(i)%display_format, i=1, num_gauges)
142 read(unit, *) (gauges(i)%min_time_increment, i=1, num_gauges)
148 allocate(gauges(i)%q_out_vars(num_eqn))
149 read(unit, *) gauges(i)%q_out_vars
152 gauges(i)%num_out_vars = 0
153 do n = 1,
size(gauges(i)%q_out_vars, 1)
154 if (gauges(i)%q_out_vars(n))
then
155 gauges(i)%num_out_vars = gauges(i)%num_out_vars + 1
161 if (num_aux > 0)
then
165 allocate(gauges(i)%aux_out_vars(num_aux))
166 read(unit, *) gauges(i)%aux_out_vars
169 do n = 1,
size(gauges(i)%aux_out_vars, 1)
170 if (gauges(i)%aux_out_vars(n))
then
171 gauges(i)%num_out_vars = gauges(i)%num_out_vars + 1
182 allocate(gauges(i)%data(gauges(i)%num_out_vars + 1, max_buffer))
187 gauges(i)%file_name =
'gaugexxxxx.txt'
188 num = gauges(i)%gauge_num
191 gauges(i)%file_name(pos:pos) = char(ichar(
'0') + digit)
197 open(unit=outgaugeunit, file=gauges(i)%file_name, &
198 status=
'old', position=
'append', form=
'formatted')
200 open(unit=outgaugeunit, file=gauges(i)%file_name, &
201 status=
'unknown', position=
'append', form=
'formatted')
205 header_1 =
"('# gauge_id= ',i5,' " // &
206 "location=( ',1e15.7,' ',1e15.7,' ) " // &
208 write(outgaugeunit, header_1) gauges(i)%gauge_num, &
211 gauges(i)%num_out_vars
216 do n=1,
size(gauges(i)%q_out_vars, 1)
217 if (gauges(i)%q_out_vars(n))
then
218 write(q_column(3 * index + 2:4 + 3 * index),
"(i3)") n
222 q_column(3 * index + 2:4 + 3 * index) =
"]"
226 if (
allocated(gauges(i)%aux_out_vars))
then
227 do n=1,
size(gauges(i)%aux_out_vars, 1)
228 if (gauges(i)%aux_out_vars(n))
then
229 write(aux_column(3 * index + 2:4 + 3 * index),
"(i3)") n
234 aux_column(3 * index + 2:4 + 3 * index) =
"]"
236 write(outgaugeunit,
"(a,a,a,a)")
"# level, time, q", &
237 trim(q_column),
", aux", &
245 module_setup = .true.
266 integer :: lev, mptr, i, k1, ki
280 if ((gauges(i)%x >= rnode(cornxlo,mptr)) .and. &
281 (gauges(i)%x <= rnode(cornxhi,mptr)) .and. &
282 (gauges(i)%y >= rnode(cornylo,mptr)) .and. &
283 (gauges(i)%y <= rnode(cornyhi,mptr)) )
then
287 mptr = node(levelptr, mptr)
294 if (mbestsrc(i) .eq. 0) &
295 print *,
"ERROR in setting grid src for gauge data", i
299 call
qsorti(mbestorder,num_gauges,mbestsrc)
319 ki = mbestsrc(mbestorder(i))
332 if (num_gauges > 0)
then
334 mbestg2(ki) = num_gauges
342 subroutine update_gauges(q, aux, xlow, ylow, num_eqn, mitot, mjtot, num_aux, &
362 use amr_module, only: nestlevel, nghost, timemult, rnode, node, maxvar
368 integer,
intent(in) :: num_eqn, mitot, mjtot, num_aux, mptr
369 real(kind=8),
intent(in) :: q(num_eqn, mitot, mjtot)
370 real(kind=8),
intent(in) :: aux(num_aux, mitot, mjtot)
371 real(kind=8),
intent(in) :: xlow, ylow
374 real(kind=8) :: var(maxvar + maxaux)
375 real(kind=8) :: xcent, ycent, xoff, yoff, tgrid, hx, hy
376 integer :: i, j, i1, i2, iindex, jindex, n, ii, index, level, var_index
379 if (num_gauges == 0)
then
392 tgrid = rnode(timemult, mptr)
393 level = node(nestlevel, mptr)
400 if (mptr /= mbestsrc(ii))
then
401 print *,
'*** should not happen... i, ii, mbestsrc(ii), mptr:'
402 print *, i, ii, mbestsrc(ii), mptr
405 if (tgrid < gauges(ii)%t_start .or. tgrid > gauges(ii)%t_end)
then
410 if (tgrid - gauges(ii)%last_time < gauges(ii)%min_time_increment)
then
416 iindex = int(.5d0 + (gauges(ii)%x - xlow) / hx)
417 jindex = int(.5d0 + (gauges(ii)%y - ylow) / hy)
418 if ((iindex < nghost .or. iindex > mitot-nghost) .or. &
419 (jindex < nghost .or. jindex > mjtot-nghost))
then
420 print *,
"ERROR in output of Gauge Data "
422 xcent = xlow + (iindex - 0.5d0) * hx
423 ycent = ylow + (jindex - 0.5d0) * hy
424 xoff = (gauges(ii)%x - xcent) / hx
425 yoff = (gauges(ii)%y - ycent) / hy
441 do n = 1,
size(gauges(ii)%q_out_vars, 1)
442 if (gauges(ii)%q_out_vars(n))
then
443 var_index = var_index + 1
444 var(var_index) = (1.d0 - xoff) * (1.d0 - yoff) * q(n, iindex, jindex) &
445 + xoff * (1.d0 - yoff) * q(n, iindex + 1, jindex) &
446 + (1.d0 - xoff) * yoff * q(n, iindex, jindex + 1) &
447 + xoff * yoff * q(n, iindex + 1, jindex + 1)
451 if (
allocated(gauges(ii)%aux_out_vars))
then
452 do n = 1,
size(gauges(ii)%aux_out_vars, 1)
453 if (gauges(ii)%aux_out_vars(n))
then
454 var_index = var_index + 1
455 var(var_index) = (1.d0 - xoff) * (1.d0 - yoff) * aux(n, iindex, jindex) &
456 + xoff * (1.d0 - yoff) * aux(n, iindex + 1, jindex) &
457 + (1.d0 - xoff) * yoff * aux(n, iindex, jindex + 1) &
458 + xoff * yoff * aux(n, iindex + 1, jindex + 1)
464 if (gauges(ii)%num_out_vars /= var_index)
then
465 print *, gauges(ii)%num_out_vars, var_index
466 print *, gauges(ii)%q_out_vars
467 print *, gauges(ii)%aux_out_vars
468 stop
"Somehow we did not grab all the values we wanted..."
472 do j = 1, gauges(ii)%num_out_vars
473 if (abs(var(j)) < 1d-90) var(j) = 0.d0
477 index = gauges(ii)%buffer_index
479 gauges(ii)%level(index) = level
480 gauges(ii)%data(1,index) = tgrid
481 do j = 1, gauges(ii)%num_out_vars
482 gauges(ii)%data(1 + j, index) = var(j)
485 gauges(ii)%buffer_index = index + 1
486 if (gauges(ii)%buffer_index > max_buffer)
then
490 gauges(ii)%last_time = tgrid
504 integer,
intent(in) :: gauge_num
507 integer :: j, k, myunit
508 integer :: omp_get_thread_num, mythread
509 character(len=32) :: out_format
514 myunit = outgaugeunit + mythread
517 if (gauges(gauge_num)%file_format == 1)
then
520 write(out_format,
"(A7, i2, A6, A1)")
"(i5.2,", &
521 gauges(gauge_num)%num_out_vars + 1, gauges(gauge_num)%display_format,
")"
523 open(unit=myunit, file=gauges(gauge_num)%file_name, status=
'old', &
524 position=
'append', form=
'formatted')
529 do j = 1, gauges(gauge_num)%buffer_index - 1
530 write(myunit, out_format) gauges(gauge_num)%level(j), &
531 (gauges(gauge_num)%data(k, j), k=1, gauges(gauge_num)%num_out_vars + 1)
533 gauges(gauge_num)%buffer_index = 1
538 print *,
"Unhandled file format ", gauges(gauge_num)%file_format
subroutine qsorti(ORD, N, A)
subroutine update_gauges(q, aux, xlow, ylow, num_eqn, mitot, mjtot, num_aux, mptr)
subroutine opendatafile(iunit, fname)
subroutine set_gauges(restart, num_eqn, num_aux, fname)
subroutine print_gauges_and_reset_nextloc(gauge_num)