4 #include "cpp_macros.h"
14 integer,
intent(in) :: n_values, n_vars
15 real(dp),
intent(inout) :: u(n_values, n_vars)
20 integer,
intent(in) :: nf
21 integer,
intent(in) :: n_var
22 integer,
intent(in) :: flux_dim
23 real(dp),
intent(in) :: u(nf, n_var)
24 real(dp),
intent(out) :: w(nf)
27 subroutine subr_flux(nf, n_var, flux_dim, u, flux, box, line_ix, s_deriv)
29 integer,
intent(in) :: nf
30 integer,
intent(in) :: n_var
31 integer,
intent(in) :: flux_dim
32 real(dp),
intent(in) :: u(nf, n_var)
33 real(dp),
intent(out) :: flux(nf, n_var)
34 type(
box_t),
intent(in) :: box
35 integer,
intent(in) :: line_ix(NDIM-1)
36 integer,
intent(in) :: s_deriv
40 n_other_dt, other_dt, box, line_ix, s_deriv)
42 integer,
intent(in) :: nf
43 integer,
intent(in) :: n_var
44 integer,
intent(in) :: flux_dim
45 real(dp),
intent(in) :: u(nf, n_var)
46 real(dp),
intent(out) :: flux(nf, n_var)
48 real(dp),
intent(out) :: cfl_sum(nf-1)
49 integer,
intent(in) :: n_other_dt
50 real(dp),
intent(inout) :: other_dt(n_other_dt)
51 type(
box_t),
intent(in) :: box
52 integer,
intent(in) :: line_ix(NDIM-1)
53 integer,
intent(in) :: s_deriv
58 integer,
intent(in) :: nf
59 integer,
intent(in) :: n_var
60 integer,
intent(in) :: flux_dim
61 real(dp),
intent(inout) :: flux(nf, n_var)
62 type(
box_t),
intent(in) :: box
63 integer,
intent(in) :: line_ix(NDIM-1)
64 integer,
intent(in) :: s_deriv
67 subroutine subr_flux_dir(box, line_ix, s_deriv, n_var, flux_dim, direction_positive)
69 type(
box_t),
intent(in) :: box
70 integer,
intent(in) :: line_ix(NDIM-1)
71 integer,
intent(in) :: s_deriv
72 integer,
intent(in) :: n_var
73 integer,
intent(in) :: flux_dim
75 logical,
intent(out) :: direction_positive(box%n_cell+1, n_var)
80 integer,
intent(in) :: n_cc
81 integer,
intent(in) :: n_var
82 real(dp),
intent(inout) :: cc_line(n_cc, n_var)
83 integer,
intent(in) :: flux_dim
84 type(
box_t),
intent(in) :: box
85 integer,
intent(in) :: line_ix(NDIM-1)
86 integer,
intent(in) :: s_deriv
89 subroutine subr_source(box, dt, n_vars, i_cc, s_deriv, s_out, n_dt, dt_lim, mask)
91 type(
box_t),
intent(inout) :: box
92 real(dp),
intent(in) :: dt
93 integer,
intent(in) :: n_vars
94 integer,
intent(in) :: i_cc(n_vars)
95 integer,
intent(in) :: s_deriv
96 integer,
intent(in) :: s_out
98 integer,
intent(in) :: n_dt
100 real(dp),
intent(inout) :: dt_lim(n_dt)
102 logical,
intent(in) :: mask(DTIMES(box%n_cell))
107 type(
box_t),
intent(in) :: box
108 logical,
intent(out) :: mask(DTIMES(box%n_cell))
112 public :: flux_diff_1d, flux_diff_2d, flux_diff_3d
113 public :: flux_koren_1d, flux_koren_2d, flux_koren_3d
114 public :: flux_upwind_1d, flux_upwind_2d, flux_upwind_3d
115 public :: flux_kurganovtadmor_1d, reconstruct_lr_1d
116 public :: flux_generic_box, flux_generic_tree
117 public :: flux_upwind_box, flux_upwind_tree
118 public :: flux_update_densities
119 public :: flux_dummy_conversion
120 public :: flux_dummy_source
121 public :: flux_dummy_modify
122 public :: flux_dummy_line_modify
123 public :: flux_get_line_cc, flux_get_line_1cc
124 public :: flux_get_line_fc, flux_get_line_1fc
129 subroutine flux_diff_1d(cc, dc, inv_dr, nc, ngc)
130 integer,
intent(in) :: nc
131 integer,
intent(in) :: ngc
132 real(dp),
intent(in) :: cc(1-ngc:nc+ngc)
134 real(dp),
intent(inout) :: dc(1:nc+1)
135 real(dp),
intent(in) :: inv_dr
139 dc(i) = -dc(i) * (cc(i) - cc(i-1)) * inv_dr
141 end subroutine flux_diff_1d
144 subroutine flux_diff_2d(cc, dc, inv_dr, nc, ngc)
145 integer,
intent(in) :: nc
146 integer,
intent(in) :: ngc
148 real(dp),
intent(in) :: cc(1-ngc:nc+ngc, 1-ngc:nc+ngc)
150 real(dp),
intent(inout) :: dc(1:nc+1, 1:nc+1, 2)
151 real(dp),
intent(in) :: inv_dr(2)
152 real(dp) :: cc_1d(1-ngc:nc+ngc), dc_1d(1:nc+1)
157 call flux_diff_1d(cc(:, n), dc(:, n, 1), inv_dr(1), nc, ngc)
162 call flux_diff_1d(cc_1d, dc_1d, inv_dr(2), nc, ngc)
165 end subroutine flux_diff_2d
168 subroutine flux_diff_3d(cc, dc, inv_dr, nc, ngc)
170 integer,
intent(in) :: nc
172 integer,
intent(in) :: ngc
174 real(dp),
intent(in) :: cc(1-ngc:nc+ngc, 1-ngc:nc+ngc, 1-ngc:nc+ngc)
176 real(dp),
intent(inout) :: dc(1:nc+1, 1:nc+1, 1:nc+1, 3)
178 real(dp),
intent(in) :: inv_dr(3)
179 real(dp) :: cc_1d(1-ngc:nc+ngc), dc_1d(1:nc+1)
185 call flux_diff_1d(cc(:, n, m), dc(:, n, m, 1), &
190 dc_1d = dc(n, :, m, 2)
191 call flux_diff_1d(cc_1d, dc_1d, inv_dr(2), nc, ngc)
192 dc(n, :, m, 2) = dc_1d
196 dc_1d = dc(n, m, :, 3)
197 call flux_diff_1d(cc_1d, dc_1d, inv_dr(3), nc, ngc)
198 dc(n, m, :, 3) = dc_1d
201 end subroutine flux_diff_3d
204 subroutine flux_koren_1d(cc, v, nc, ngc)
205 integer,
intent(in) :: nc
206 integer,
intent(in) :: ngc
207 real(dp),
intent(in) :: cc(1-ngc:nc+ngc)
209 real(dp),
intent(inout) :: v(1:nc+1)
210 real(dp) :: gradp, gradc, gradn
214 gradc = cc(n) - cc(n-1)
215 if (v(n) < 0.0_dp)
then
216 gradn = cc(n+1) - cc(n)
217 v(n) = v(n) * (cc(n) - 0.5_dp * af_limiter_koren(gradc, gradn))
219 gradp = cc(n-1) - cc(n-2)
220 v(n) = v(n) * (cc(n-1) + 0.5_dp * af_limiter_koren(gradc, gradp))
224 end subroutine flux_koren_1d
227 subroutine flux_koren_2d(cc, v, nc, ngc)
228 integer,
intent(in) :: nc
229 integer,
intent(in) :: ngc
231 real(dp),
intent(in) :: cc(1-ngc:nc+ngc, 1-ngc:nc+ngc)
233 real(dp),
intent(inout) :: v(1:nc+1, 1:nc+1, 2)
234 real(dp) :: cc_1d(1-ngc:nc+ngc), v_1d(1:nc+1)
239 call flux_koren_1d(cc(:, n), v(:, n, 1), nc, ngc)
244 call flux_koren_1d(cc_1d, v_1d, nc, ngc)
247 end subroutine flux_koren_2d
252 subroutine reconstruct_lr_1d(nc, ngc, n_vars, cc, u_l, u_r, limiter)
253 integer,
intent(in) :: nc
254 integer,
intent(in) :: ngc
255 integer,
intent(in) :: n_vars
256 real(dp),
intent(in) :: cc(1-ngc:nc+ngc, n_vars)
258 real(dp),
intent(inout) :: u_l(1:nc+1, n_vars)
260 real(dp),
intent(inout) :: u_r(1:nc+1, n_vars)
261 integer,
intent(in) :: limiter
262 real(dp) :: slopes(0:nc+1, n_vars)
264 associate(a=>cc(1:nc+2, :) - cc(0:nc+1, :), &
265 b=>cc(0:nc+1, :) - cc(-1:nc, :))
267 slopes = af_limiter_apply(a, b, limiter)
270 u_l(1:nc+1, :) = cc(0:nc, :) + 0.5_dp * slopes(0:nc, :)
272 if (af_limiter_symmetric(limiter))
then
273 u_r(1:nc+1, :) = cc(1:nc+1, :) - 0.5_dp * slopes(1:nc+1, :)
275 slopes = af_limiter_apply(b, a, limiter)
276 u_r(1:nc+1, :) = cc(1:nc+1, :) - 0.5_dp * slopes(1:nc+1, :)
279 end subroutine reconstruct_lr_1d
282 subroutine reconstruct_upwind_1d(nc, ngc, n_vars, cc, u_f, limiter, direction_positive)
283 integer,
intent(in) :: nc
284 integer,
intent(in) :: ngc
285 integer,
intent(in) :: n_vars
286 real(dp),
intent(in) :: cc(1-ngc:nc+ngc, n_vars)
288 real(dp),
intent(inout) :: u_f(1:nc+1, n_vars)
289 integer,
intent(in) :: limiter
291 logical,
intent(in) :: direction_positive(nc+1, n_vars)
293 associate(a=>cc(1:nc+2, :) - cc(0:nc+1, :), &
294 b=>cc(0:nc+1, :) - cc(-1:nc, :))
295 where (direction_positive)
296 u_f = cc(0:nc, :) + 0.5_dp * &
297 af_limiter_apply(a(1:nc+1, :), b(1:nc+1, :), limiter)
299 u_f = cc(1:nc+1, :) - 0.5_dp * &
300 af_limiter_apply(b(2:nc+2, :), a(2:nc+2, :), limiter)
303 end subroutine reconstruct_upwind_1d
306 subroutine flux_kurganovtadmor_1d(n_values, n_vars, flux_l, flux_r, &
307 u_l, u_r, wmax, flux)
308 integer,
intent(in) :: n_values
309 integer,
intent(in) :: n_vars
310 real(dp),
intent(in) :: flux_l(n_values, n_vars)
311 real(dp),
intent(in) :: flux_r(n_values, n_vars)
312 real(dp),
intent(in) :: u_l(n_values, n_vars)
313 real(dp),
intent(in) :: u_r(n_values, n_vars)
314 real(dp),
intent(in) :: wmax(n_values)
315 real(dp),
intent(inout) :: flux(n_values, n_vars)
317 flux = 0.5_dp * (flux_l + flux_r - spread(wmax, 2, n_vars) * (u_r - u_l))
318 end subroutine flux_kurganovtadmor_1d
320 subroutine flux_update_densities(tree, dt, n_vars, i_cc, n_vars_flux, i_cc_flux, i_flux, &
321 s_deriv, n_prev, s_prev, w_prev, s_out, add_source_box, n_dt, dt_lim, get_mask)
322 type(af_t),
intent(inout) :: tree
323 real(dp),
intent(in) :: dt
324 integer,
intent(in) :: n_vars
325 integer,
intent(in) :: i_cc(n_vars)
326 integer,
intent(in) :: n_vars_flux
327 integer,
intent(in) :: i_cc_flux(n_vars_flux)
328 integer,
intent(in) :: i_flux(n_vars_flux)
329 integer,
intent(in) :: s_deriv
330 integer,
intent(in) :: n_prev
331 integer,
intent(in) :: s_prev(n_prev)
332 real(dp),
intent(in) :: w_prev(n_prev)
333 integer,
intent(in) :: s_out
337 integer,
intent(in) :: n_dt
339 real(dp),
intent(out) :: dt_lim(n_dt)
342 integer :: lvl, n, m, id, ijk, nc, i_var, iv
343 logical :: mask(dtimes(tree%n_cell))
344 real(dp) :: tmp(dtimes(tree%n_cell))
345 real(dp) :: dt_dr(ndim), my_dt(n_dt)
346 real(dp) :: rfac(2, tree%n_cell)
357 do lvl = 1, tree%highest_lvl
359 do n = 1,
size(tree%lvls(lvl)%leaves)
360 id = tree%lvls(lvl)%leaves(n)
361 dt_dr = dt/tree%boxes(id)%dr
363 associate(cc => tree%boxes(id)%cc, fc => tree%boxes(id)%fc)
364 if (
present(get_mask))
then
365 call get_mask(tree%boxes(id), mask)
376 tmp = tmp + w_prev(m) * cc(dtimes(1:nc), iv+s_prev(m))
379 cc(dtimes(1:nc), iv+s_out) = tmp
382 call add_source_box(tree%boxes(id), dt, n_vars, i_cc, s_deriv, &
383 s_out, n_dt, my_dt, mask)
384 dt_lim = min(dt_lim, my_dt)
389 cc(i, i_cc_flux+s_out) = cc(i, i_cc_flux+s_out) + &
391 (fc(i, 1, i_flux) - fc(i+1, 1, i_flux))
395 if (tree%coord_t == af_cyl)
then
396 call af_cyl_flux_factors(tree%boxes(id), rfac)
399 cc(i, j, i_cc_flux+s_out) = cc(i, j, i_cc_flux+s_out) + &
401 rfac(1, i) * fc(i, j, 1, i_flux) - &
402 rfac(2, i) * fc(i+1, j, 1, i_flux)) &
404 (fc(i, j, 2, i_flux) - fc(i, j+1, 2, i_flux))
410 cc(i, j, i_cc_flux+s_out) = cc(i, j, i_cc_flux+s_out) + &
412 (fc(i, j, 1, i_flux) - fc(i+1, j, 1, i_flux)) &
414 (fc(i, j, 2, i_flux) - fc(i, j+1, 2, i_flux))
421 cc(i, j, k, i_cc_flux+s_out) = cc(i, j, k, i_cc_flux+s_out) + &
423 (fc(i, j, k, 1, i_flux) - fc(i+1, j, k, 1, i_flux)) + &
425 (fc(i, j, k, 2, i_flux) - fc(i, j+1, k, 2, i_flux)) + &
427 (fc(i, j, k, 3, i_flux) - fc(i, j, k+1, 3, i_flux))
436 end subroutine flux_update_densities
439 subroutine flux_generic_tree(tree, n_vars, i_cc, s_deriv, i_flux, dt_lim, &
440 max_wavespeed, flux_from_primitives, flux_modify, line_modify, &
441 to_primitive, to_conservative, limiter)
444 type(af_t),
intent(inout) :: tree
445 integer,
intent(in) :: n_vars
446 integer,
intent(in) :: i_cc(n_vars)
447 integer,
intent(in) :: s_deriv
448 integer,
intent(in) :: i_flux(n_vars)
450 real(dp),
intent(out) :: dt_lim
454 procedure(
subr_flux) :: flux_from_primitives
464 integer,
intent(in) :: limiter
470 call af_restrict_ref_boundary(tree, i_cc+s_deriv)
476 do lvl = 1, tree%highest_lvl
478 do i = 1,
size(tree%lvls(lvl)%leaves)
479 call flux_generic_box(tree, tree%lvls(lvl)%leaves(i), tree%n_cell, &
480 n_vars, i_cc, s_deriv, i_flux, my_dt, max_wavespeed, &
481 flux_from_primitives, flux_modify, line_modify, &
482 to_primitive, to_conservative, limiter)
483 dt_lim = min(dt_lim, my_dt)
490 call af_consistent_fluxes(tree, i_flux)
492 end subroutine flux_generic_tree
495 subroutine flux_generic_box(tree, id, nc, n_vars, i_cc, s_deriv, i_flux, dt_lim, &
496 max_wavespeed, flux_from_primitives, flux_modify, line_modify, &
497 to_primitive, to_conservative, limiter)
500 type(
af_t),
intent(inout) :: tree
501 integer,
intent(in) :: id
502 integer,
intent(in) :: nc
503 integer,
intent(in) :: n_vars
504 integer,
intent(in) :: i_cc(n_vars)
505 integer,
intent(in) :: s_deriv
506 integer,
intent(in) :: i_flux(n_vars)
508 real(dp),
intent(out) :: dt_lim
512 procedure(
subr_flux) :: flux_from_primitives
522 integer,
intent(in) :: limiter
525 real(dp) :: cc(dtimes(-1:nc+2), n_vars)
526 real(dp) :: cc_line(-1:nc+2, n_vars)
527 real(dp) :: flux(nc+1, n_vars)
528 real(dp) :: u_l(nc+1, n_vars), u_r(nc+1, n_vars)
529 real(dp) :: w_l(nc+1), w_r(nc+1)
530 real(dp) :: flux_l(nc+1, n_vars), flux_r(nc+1, n_vars)
531 real(dp) :: cfl_sum(dtimes(nc)), cfl_sum_line(nc), inv_dr(ndim)
532 integer :: flux_dim, line_ix(ndim-1)
540 inv_dr = 1/tree%boxes(id)%dr
543 call af_gc2_box(tree, id, i_cc+s_deriv, cc)
554 do flux_dim = 1, ndim
562 select case (flux_dim)
568 cc_line = cc(:, i, :)
570 cc_line = cc(i, :, :)
573 cc_line = cc(:, i, j, :)
575 cc_line = cc(i, :, j, :)
577 cc_line = cc(i, j, :, :)
588 call line_modify(nc+4, n_vars, cc_line, flux_dim, &
589 tree%boxes(id), line_ix, s_deriv)
593 call to_primitive(nc+4, n_vars, cc_line)
597 call reconstruct_lr_1d(nc, 2, n_vars, cc_line, u_l, u_r, limiter)
600 call max_wavespeed(nc+1, n_vars, flux_dim, u_l, w_l)
601 call max_wavespeed(nc+1, n_vars, flux_dim, u_r, w_r)
604 call flux_from_primitives(nc+1, n_vars, flux_dim, u_l, flux_l, &
605 tree%boxes(id), line_ix, s_deriv)
606 call flux_from_primitives(nc+1, n_vars, flux_dim, u_r, flux_r, &
607 tree%boxes(id), line_ix, s_deriv)
610 call to_conservative(nc+1, n_vars, u_l)
611 call to_conservative(nc+1, n_vars, u_r)
617 cfl_sum_line = max(w_l(1:nc), w_l(2:nc+1)) * inv_dr(ndim)
620 call flux_kurganovtadmor_1d(nc+1, n_vars, flux_l, flux_r, &
624 call flux_modify(nc+1, n_vars, flux_dim, flux, &
625 tree%boxes(id), line_ix, s_deriv)
628 select case (flux_dim)
631 tree%boxes(id)%fc(:, flux_dim, i_flux) = flux
632 cfl_sum = cfl_sum + cfl_sum_line
635 tree%boxes(id)%fc(:, i, flux_dim, i_flux) = flux
636 cfl_sum(:, i) = cfl_sum(:, i) + cfl_sum_line
638 tree%boxes(id)%fc(i, :, flux_dim, i_flux) = flux
639 cfl_sum(i, :) = cfl_sum(i, :) + cfl_sum_line
642 tree%boxes(id)%fc(:, i, j, flux_dim, i_flux) = flux
643 cfl_sum(:, i, j) = cfl_sum(:, i, j) + cfl_sum_line
645 tree%boxes(id)%fc(i, :, j, flux_dim, i_flux) = flux
646 cfl_sum(i, :, j) = cfl_sum(i, :, j) + cfl_sum_line
648 tree%boxes(id)%fc(i, j, :, flux_dim, i_flux) = flux
649 cfl_sum(i, j, :) = cfl_sum(i, j, :) + cfl_sum_line
661 dt_lim = 1/max(maxval(cfl_sum), 1e-100_dp)
663 end subroutine flux_generic_box
666 subroutine flux_upwind_tree(tree, n_vars, i_cc, s_deriv, i_flux, n_dt, dt_lim, &
667 flux_upwind, flux_direction, line_modify, limiter)
670 type(af_t),
intent(inout) :: tree
671 integer,
intent(in) :: n_vars
672 integer,
intent(in) :: i_cc(n_vars)
673 integer,
intent(in) :: s_deriv
674 integer,
intent(in) :: i_flux(n_vars)
676 integer,
intent(in) :: n_dt
678 real(dp),
intent(out) :: dt_lim(n_dt)
686 integer,
intent(in) :: limiter
688 real(dp) :: my_dt(n_dt)
691 call af_restrict_ref_boundary(tree, i_cc+s_deriv)
697 do lvl = 1, tree%highest_lvl
699 do i = 1,
size(tree%lvls(lvl)%leaves)
700 call flux_upwind_box(tree, tree%lvls(lvl)%leaves(i), tree%n_cell, &
701 n_vars, i_cc, s_deriv, i_flux, n_dt, my_dt, flux_upwind, &
702 flux_direction, line_modify, limiter)
703 dt_lim = min(dt_lim, my_dt)
710 call af_consistent_fluxes(tree, i_flux)
712 end subroutine flux_upwind_tree
715 subroutine flux_upwind_box(tree, id, nc, n_vars, i_cc, s_deriv, i_flux, &
716 n_dt, dt_lim, flux_upwind, flux_direction, line_modify, limiter)
719 type(
af_t),
intent(inout) :: tree
720 integer,
intent(in) :: id
721 integer,
intent(in) :: nc
722 integer,
intent(in) :: n_vars
723 integer,
intent(in) :: i_cc(n_vars)
724 integer,
intent(in) :: s_deriv
725 integer,
intent(in) :: i_flux(n_vars)
727 integer,
intent(in) :: n_dt
728 real(dp),
intent(inout) :: dt_lim(n_dt)
736 integer,
intent(in) :: limiter
738 real(dp) :: cc(dtimes(-1:nc+2), n_vars)
739 real(dp) :: cc_line(-1:nc+2, n_vars)
740 real(dp) :: flux(nc+1, n_vars)
741 real(dp) :: u_l(nc+1, n_vars)
742 real(dp) :: cfl_sum(dtimes(nc)), cfl_sum_line(nc)
743 real(dp) :: other_dt(n_dt-1)
744 logical :: direction_positive(nc+1, n_vars)
745 integer :: flux_dim, line_ix(ndim-1)
754 call af_gc2_box(tree, id, i_cc+s_deriv, cc)
760 dt_lim(2:) = 1e100_dp
763 associate(fc => tree%boxes(id)%fc)
764 do flux_dim = 1, ndim
772 select case (flux_dim)
778 cc_line = cc(:, i, :)
780 cc_line = cc(i, :, :)
783 cc_line = cc(:, i, j, :)
785 cc_line = cc(i, :, j, :)
787 cc_line = cc(i, j, :, :)
796 call flux_direction(tree%boxes(id), line_ix, s_deriv, &
797 n_vars, flux_dim, direction_positive)
800 call line_modify(nc+4, n_vars, cc_line, flux_dim, &
801 tree%boxes(id), line_ix, s_deriv)
804 call reconstruct_upwind_1d(nc, 2, n_vars, cc_line, u_l, &
805 limiter, direction_positive)
807 call flux_upwind(nc+1, n_vars, flux_dim, u_l, flux, cfl_sum_line, &
808 n_dt-1, other_dt, tree%boxes(id), line_ix, s_deriv)
809 dt_lim(2:) = min(dt_lim(2:), other_dt)
812 select case (flux_dim)
815 fc(:, flux_dim, i_flux) = flux
816 cfl_sum = cfl_sum + cfl_sum_line
819 fc(:, i, flux_dim, i_flux) = flux
820 cfl_sum(:, i) = cfl_sum(:, i) + cfl_sum_line
822 fc(i, :, flux_dim, i_flux) = flux
823 cfl_sum(i, :) = cfl_sum(i, :) + cfl_sum_line
826 fc(:, i, j, flux_dim, i_flux) = flux
827 cfl_sum(:, i, j) = cfl_sum(:, i, j) + cfl_sum_line
829 fc(i, :, j, flux_dim, i_flux) = flux
830 cfl_sum(i, :, j) = cfl_sum(i, :, j) + cfl_sum_line
832 fc(i, j, :, flux_dim, i_flux) = flux
833 cfl_sum(i, j, :) = cfl_sum(i, j, :) + cfl_sum_line
846 dt_lim(1) = 1/maxval(cfl_sum)
848 end subroutine flux_upwind_box
853 subroutine flux_get_line_cc(box, ivs, flux_dim, line_ix, cc_line)
854 type(box_t),
intent(in) :: box
855 integer,
intent(in) :: ivs(:)
856 integer,
intent(in) :: flux_dim
857 integer,
intent(in) :: line_ix(ndim-1)
858 real(dp),
intent(inout) :: cc_line(box%n_cell+2, size(ivs))
860 select case (flux_dim)
863 cc_line = box%cc(:, ivs)
866 cc_line = box%cc(:, line_ix(1), ivs)
868 cc_line = box%cc(line_ix(1), :, ivs)
871 cc_line = box%cc(:, line_ix(1), line_ix(2), ivs)
873 cc_line = box%cc(line_ix(1), :, line_ix(2), ivs)
875 cc_line = box%cc(line_ix(1), line_ix(2), :, ivs)
878 end subroutine flux_get_line_cc
883 subroutine flux_get_line_1cc(box, iv, flux_dim, line_ix, cc_line)
884 type(box_t),
intent(in) :: box
885 integer,
intent(in) :: iv
886 integer,
intent(in) :: flux_dim
887 integer,
intent(in) :: line_ix(ndim-1)
888 real(dp),
intent(inout) :: cc_line(box%n_cell+2)
890 select case (flux_dim)
893 cc_line = box%cc(:, iv)
896 cc_line = box%cc(:, line_ix(1), iv)
898 cc_line = box%cc(line_ix(1), :, iv)
901 cc_line = box%cc(:, line_ix(1), line_ix(2), iv)
903 cc_line = box%cc(line_ix(1), :, line_ix(2), iv)
905 cc_line = box%cc(line_ix(1), line_ix(2), :, iv)
908 end subroutine flux_get_line_1cc
912 subroutine flux_get_line_fc(box, ivs, flux_dim, line_ix, fc_line)
913 type(box_t),
intent(in) :: box
914 integer,
intent(in) :: ivs(:)
915 integer,
intent(in) :: flux_dim
916 integer,
intent(in) :: line_ix(ndim-1)
917 real(dp),
intent(inout) :: fc_line(box%n_cell+1, size(ivs))
919 select case (flux_dim)
922 fc_line = box%fc(:, flux_dim, ivs)
925 fc_line = box%fc(:, line_ix(1), flux_dim, ivs)
927 fc_line = box%fc(line_ix(1), :, flux_dim, ivs)
930 fc_line = box%fc(:, line_ix(1), line_ix(2), flux_dim, ivs)
932 fc_line = box%fc(line_ix(1), :, line_ix(2), flux_dim, ivs)
934 fc_line = box%fc(line_ix(1), line_ix(2), :, flux_dim, ivs)
937 end subroutine flux_get_line_fc
941 subroutine flux_get_line_1fc(box, iv, flux_dim, line_ix, fc_line)
942 type(box_t),
intent(in) :: box
943 integer,
intent(in) :: iv
944 integer,
intent(in) :: flux_dim
945 integer,
intent(in) :: line_ix(ndim-1)
946 real(dp),
intent(inout) :: fc_line(box%n_cell+1)
948 select case (flux_dim)
951 fc_line = box%fc(:, flux_dim, iv)
954 fc_line = box%fc(:, line_ix(1), flux_dim, iv)
956 fc_line = box%fc(line_ix(1), :, flux_dim, iv)
959 fc_line = box%fc(:, line_ix(1), line_ix(2), flux_dim, iv)
961 fc_line = box%fc(line_ix(1), :, line_ix(2), flux_dim, iv)
963 fc_line = box%fc(line_ix(1), line_ix(2), :, flux_dim, iv)
966 end subroutine flux_get_line_1fc
969 subroutine flux_koren_3d(cc, v, nc, ngc)
971 integer,
intent(in) :: nc
973 integer,
intent(in) :: ngc
975 real(dp),
intent(in) :: cc(1-ngc:nc+ngc, 1-ngc:nc+ngc, 1-ngc:nc+ngc)
977 real(dp),
intent(inout) :: v(1:nc+1, 1:nc+1, 1:nc+1, 3)
978 real(dp) :: cc_1d(1-ngc:nc+ngc), v_1d(1:nc+1)
984 call flux_koren_1d(cc(:, n, m), v(:, n, m, 1), &
990 call flux_koren_1d(cc_1d, v_1d, nc, ngc)
996 call flux_koren_1d(cc_1d, v_1d, nc, ngc)
1000 end subroutine flux_koren_3d
1003 subroutine flux_upwind_1d(cc, v, nc, ngc)
1004 integer,
intent(in) :: nc
1005 integer,
intent(in) :: ngc
1006 real(dp),
intent(in) :: cc(1-ngc:nc+ngc)
1008 real(dp),
intent(inout) :: v(1:nc+1)
1012 if (v(n) < 0.0_dp)
then
1015 v(n) = v(n) * cc(n-1)
1019 end subroutine flux_upwind_1d
1022 subroutine flux_upwind_2d(cc, v, nc, ngc)
1023 integer,
intent(in) :: nc
1024 integer,
intent(in) :: ngc
1026 real(dp),
intent(in) :: cc(1-ngc:nc+ngc, 1-ngc:nc+ngc)
1028 real(dp),
intent(inout) :: v(1:nc+1, 1:nc+1, 2)
1029 real(dp) :: cc_1d(1-ngc:nc+ngc), v_1d(1:nc+1)
1034 call flux_upwind_1d(cc(:, n), v(:, n, 1), nc, ngc)
1039 call flux_upwind_1d(cc_1d, v_1d, nc, ngc)
1042 end subroutine flux_upwind_2d
1045 subroutine flux_upwind_3d(cc, v, nc, ngc)
1047 integer,
intent(in) :: nc
1049 integer,
intent(in) :: ngc
1051 real(dp),
intent(in) :: cc(1-ngc:nc+ngc, 1-ngc:nc+ngc, 1-ngc:nc+ngc)
1053 real(dp),
intent(inout) :: v(1:nc+1, 1:nc+1, 1:nc+1, 3)
1054 real(dp) :: cc_1d(1-ngc:nc+ngc), v_1d(1:nc+1)
1060 call flux_upwind_1d(cc(:, n, m), v(:, n, m, 1), &
1065 v_1d = v(n, :, m, 2)
1066 call flux_upwind_1d(cc_1d, v_1d, nc, ngc)
1067 v(n, :, m, 2) = v_1d
1071 v_1d = v(n, m, :, 3)
1072 call flux_upwind_1d(cc_1d, v_1d, nc, ngc)
1073 v(n, m, :, 3) = v_1d
1076 end subroutine flux_upwind_3d
1079 subroutine flux_dummy_conversion(n_values, n_vars, u)
1080 integer,
intent(in) :: n_values, n_vars
1081 real(dp),
intent(inout) :: u(n_values, n_vars)
1082 end subroutine flux_dummy_conversion
1084 subroutine flux_dummy_source(box, dt, n_vars, i_cc, s_deriv, s_out, n_dt, dt_lim, mask)
1085 type(box_t),
intent(inout) :: box
1086 real(dp),
intent(in) :: dt
1087 integer,
intent(in) :: n_vars
1088 integer,
intent(in) :: i_cc(n_vars)
1089 integer,
intent(in) :: s_deriv
1090 integer,
intent(in) :: s_out
1091 integer,
intent(in) :: n_dt
1092 real(dp),
intent(inout) :: dt_lim(n_dt)
1093 logical,
intent(in) :: mask(dtimes(box%n_cell))
1094 end subroutine flux_dummy_source
1096 subroutine flux_dummy_modify(nf, n_var, flux_dim, flux, box, line_ix, s_deriv)
1097 integer,
intent(in) :: nf
1098 integer,
intent(in) :: n_var
1099 integer,
intent(in) :: flux_dim
1100 real(dp),
intent(inout) :: flux(nf, n_var)
1101 type(box_t),
intent(in) :: box
1102 integer,
intent(in) :: line_ix(ndim-1)
1103 integer,
intent(in) :: s_deriv
1104 end subroutine flux_dummy_modify
1106 subroutine flux_dummy_line_modify(n_cc, n_var, cc_line, flux_dim, box, &
1108 integer,
intent(in) :: n_cc
1109 integer,
intent(in) :: n_var
1110 real(dp),
intent(inout) :: cc_line(n_cc, n_var)
1111 integer,
intent(in) :: flux_dim
1112 type(box_t),
intent(in) :: box
1113 integer,
intent(in) :: line_ix(ndim-1)
1114 integer,
intent(in) :: s_deriv
1115 end subroutine flux_dummy_line_modify
This module contains the core routines of Afivo, namely those that deal with initializing and changin...
Module containing a couple flux schemes for solving hyperbolic problems explicitly,...
This module contains routines related to the filling of ghost cells. Note that corner ghost cells are...
Module containing slope limiters.
This module contains routines for restriction: going from fine to coarse variables.
This module contains the basic types and constants that are used in the NDIM-dimensional version of A...
Type which stores all the boxes and levels, as well as some information about the number of boxes,...
The basic building block of afivo: a box with cell-centered and face centered data,...