Afivo 0.3
All Classes Namespaces Functions Variables Pages
m_af_flux_schemes.f90
1!> Module containing a couple flux schemes for solving hyperbolic problems
2!> explicitly, as well as handling diffusion explicitly.
4#include "cpp_macros.h"
5 use m_af_types
7
8 implicit none
9 private
10
11 interface
12 subroutine subr_prim_cons(n_values, n_vars, u)
13 import
14 integer, intent(in) :: n_values, n_vars
15 real(dp), intent(inout) :: u(n_values, n_vars)
16 end subroutine subr_prim_cons
17
18 subroutine subr_max_wavespeed(nf, n_var, flux_dim, u, w)
19 import
20 integer, intent(in) :: nf !< Number of cell faces
21 integer, intent(in) :: n_var !< Number of variables
22 integer, intent(in) :: flux_dim !< In which dimension fluxes are computed
23 real(dp), intent(in) :: u(nf, n_var)
24 real(dp), intent(out) :: w(nf)
25 end subroutine subr_max_wavespeed
26
27 subroutine subr_flux(nf, n_var, flux_dim, u, flux, box, line_ix, s_deriv)
28 import
29 integer, intent(in) :: nf !< Number of cell faces
30 integer, intent(in) :: n_var !< Number of variables
31 integer, intent(in) :: flux_dim !< In which dimension fluxes are computed
32 real(dp), intent(in) :: u(nf, n_var) !< Primitive variables
33 real(dp), intent(out) :: flux(nf, n_var) !< Computed fluxes
34 type(box_t), intent(in) :: box !< Current box
35 integer, intent(in) :: line_ix(NDIM-1) !< Index of line for dim /= flux_dim
36 integer, intent(in) :: s_deriv !< State to compute derivatives from
37 end subroutine subr_flux
38
39 subroutine subr_flux_upwind(nf, n_var, flux_dim, u, flux, cfl_sum, &
40 n_other_dt, other_dt, box, line_ix, s_deriv)
41 import
42 integer, intent(in) :: nf !< Number of cell faces
43 integer, intent(in) :: n_var !< Number of variables
44 integer, intent(in) :: flux_dim !< In which dimension fluxes are computed
45 real(dp), intent(in) :: u(nf, n_var) !< Face values
46 real(dp), intent(out) :: flux(nf, n_var) !< Computed fluxes
47 !> Terms per cell-center to be added to CFL sum, see flux_upwind_box
48 real(dp), intent(out) :: cfl_sum(nf-1)
49 integer, intent(in) :: n_other_dt !< Number of non-cfl time step restrictions
50 real(dp), intent(inout) :: other_dt(n_other_dt) !< Non-cfl time step restrictions
51 type(box_t), intent(in) :: box !< Current box
52 integer, intent(in) :: line_ix(NDIM-1) !< Index of line for dim /= flux_dim
53 integer, intent(in) :: s_deriv !< State to compute derivatives from
54 end subroutine subr_flux_upwind
55
56 subroutine subr_flux_modify(nf, n_var, flux_dim, flux, box, line_ix, s_deriv)
57 import
58 integer, intent(in) :: nf !< Number of cell faces
59 integer, intent(in) :: n_var !< Number of variables
60 integer, intent(in) :: flux_dim !< In which dimension fluxes are computed
61 real(dp), intent(inout) :: flux(nf, n_var) !< Flux that will be modified
62 type(box_t), intent(in) :: box !< Current box
63 integer, intent(in) :: line_ix(NDIM-1) !< Index of line for dim /= flux_dim
64 integer, intent(in) :: s_deriv !< State to compute derivatives from
65 end subroutine subr_flux_modify
66
67 subroutine subr_flux_dir(box, line_ix, s_deriv, n_var, flux_dim, direction_positive)
68 import
69 type(box_t), intent(in) :: box !< Current box
70 integer, intent(in) :: line_ix(NDIM-1) !< Index of line for dim /= flux_dim
71 integer, intent(in) :: s_deriv !< State to compute derivatives from
72 integer, intent(in) :: n_var !< Number of variables
73 integer, intent(in) :: flux_dim !< In which dimension fluxes are computed
74 !> True means positive flow (to the "right"), false to the left
75 logical, intent(out) :: direction_positive(box%n_cell+1, n_var)
76 end subroutine subr_flux_dir
77
78 subroutine subr_line_modify(n_cc, n_var, cc_line, flux_dim, box, line_ix, s_deriv)
79 import
80 integer, intent(in) :: n_cc !< Number of cell centers
81 integer, intent(in) :: n_var !< Number of variables
82 real(dp), intent(inout) :: cc_line(n_cc, n_var) !< Line values to modify
83 integer, intent(in) :: flux_dim !< In which dimension fluxes are computed
84 type(box_t), intent(in) :: box !< Current box
85 integer, intent(in) :: line_ix(NDIM-1) !< Index of line for dim /= flux_dim
86 integer, intent(in) :: s_deriv !< State to compute derivatives from
87 end subroutine subr_line_modify
88
89 subroutine subr_source(box, dt, n_vars, i_cc, s_deriv, s_out, n_dt, dt_lim, mask)
90 import
91 type(box_t), intent(inout) :: box
92 real(dp), intent(in) :: dt !< Time step
93 integer, intent(in) :: n_vars !< Number of cell-centered variables
94 integer, intent(in) :: i_cc(n_vars) !< Indices of all cell-centered variables
95 integer, intent(in) :: s_deriv !< State to compute derivatives from
96 integer, intent(in) :: s_out !< State to update
97 !> Number of time steps restrictions
98 integer, intent(in) :: n_dt
99 !> Maximal allowed time steps
100 real(dp), intent(inout) :: dt_lim(n_dt)
101 !> Mask where to update solution
102 logical, intent(in) :: mask(DTIMES(box%n_cell))
103 end subroutine subr_source
104
105 subroutine subr_box_mask(box, mask)
106 import
107 type(box_t), intent(in) :: box
108 logical, intent(out) :: mask(DTIMES(box%n_cell))
109 end subroutine subr_box_mask
110 end interface
111
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
125
126contains
127
128 !> Compute diffusive flux (second order)
129 subroutine flux_diff_1d(cc, dc, inv_dr, nc, ngc)
130 integer, intent(in) :: nc !< Number of cells
131 integer, intent(in) :: ngc !< Number of ghost cells
132 real(dp), intent(in) :: cc(1-ngc:nc+ngc) !< Cell-centered values
133 !> Input: diffusion coeff. at interfaces, output: fluxes
134 real(dp), intent(inout) :: dc(1:nc+1)
135 real(dp), intent(in) :: inv_dr !< Inverse grid spacing
136 integer :: i
137
138 do i = 1, nc+1
139 dc(i) = -dc(i) * (cc(i) - cc(i-1)) * inv_dr
140 end do
141 end subroutine flux_diff_1d
142
143 !> Compute diffusive flux (second order)
144 subroutine flux_diff_2d(cc, dc, inv_dr, nc, ngc)
145 integer, intent(in) :: nc !< Number of cells
146 integer, intent(in) :: ngc !< Number of ghost cells
147 !> Cell-centered values
148 real(dp), intent(in) :: cc(1-ngc:nc+ngc, 1-ngc:nc+ngc)
149 !> Input: diffusion coeff. at interfaces, output: fluxes
150 real(dp), intent(inout) :: dc(1:nc+1, 1:nc+1, 2)
151 real(dp), intent(in) :: inv_dr(2) !< Inverse grid spacing
152 real(dp) :: cc_1d(1-ngc:nc+ngc), dc_1d(1:nc+1)
153 integer :: n
154
155 do n = 1, nc
156 ! x-fluxes
157 call flux_diff_1d(cc(:, n), dc(:, n, 1), inv_dr(1), nc, ngc)
158
159 ! y-fluxes (use temporary variables for efficiency)
160 cc_1d = cc(n, :)
161 dc_1d = dc(n, :, 2)
162 call flux_diff_1d(cc_1d, dc_1d, inv_dr(2), nc, ngc)
163 dc(n, :, 2) = dc_1d ! Copy result
164 end do
165 end subroutine flux_diff_2d
166
167 !> Compute diffusive flux (second order)
168 subroutine flux_diff_3d(cc, dc, inv_dr, nc, ngc)
169 !> Number of cells
170 integer, intent(in) :: nc
171 !> Number of ghost cells
172 integer, intent(in) :: ngc
173 !> Cell-centered values
174 real(dp), intent(in) :: cc(1-ngc:nc+ngc, 1-ngc:nc+ngc, 1-ngc:nc+ngc)
175 !> Input: diffusion coeff. at interfaces, output: fluxes
176 real(dp), intent(inout) :: dc(1:nc+1, 1:nc+1, 1:nc+1, 3)
177 !> Inverse grid spacing
178 real(dp), intent(in) :: inv_dr(3)
179 real(dp) :: cc_1d(1-ngc:nc+ngc), dc_1d(1:nc+1)
180 integer :: n, m
181
182 do m = 1, nc
183 do n = 1, nc
184 ! x-fluxes
185 call flux_diff_1d(cc(:, n, m), dc(:, n, m, 1), &
186 inv_dr(1), nc, ngc)
187
188 ! y-fluxes (use temporary variables for efficiency)
189 cc_1d = cc(n, :, m)
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 ! Copy result
193
194 ! z-fluxes (use temporary variables for efficiency)
195 cc_1d = cc(n, m, :)
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 ! Copy result
199 end do
200 end do
201 end subroutine flux_diff_3d
202
203 !> Compute flux according to Koren limiter
204 subroutine flux_koren_1d(cc, v, nc, ngc)
205 integer, intent(in) :: nc !< Number of cells
206 integer, intent(in) :: ngc !< Number of ghost cells
207 real(dp), intent(in) :: cc(1-ngc:nc+ngc) !< Cell-centered values
208 !> Input: velocities at interfaces, output: fluxes
209 real(dp), intent(inout) :: v(1:nc+1)
210 real(dp) :: gradp, gradc, gradn
211 integer :: n
212
213 do n = 1, nc+1
214 gradc = cc(n) - cc(n-1) ! Current gradient
215 if (v(n) < 0.0_dp) then
216 gradn = cc(n+1) - cc(n) ! Next gradient
217 v(n) = v(n) * (cc(n) - 0.5_dp * af_limiter_koren(gradc, gradn))
218 else ! v(n) > 0
219 gradp = cc(n-1) - cc(n-2) ! Previous gradient
220 v(n) = v(n) * (cc(n-1) + 0.5_dp * af_limiter_koren(gradc, gradp))
221 end if
222 end do
223
224 end subroutine flux_koren_1d
225
226 !> Compute flux according to Koren limiter
227 subroutine flux_koren_2d(cc, v, nc, ngc)
228 integer, intent(in) :: nc !< Number of cells
229 integer, intent(in) :: ngc !< Number of ghost cells
230 !> Cell-centered values
231 real(dp), intent(in) :: cc(1-ngc:nc+ngc, 1-ngc:nc+ngc)
232 !> Input: velocities at interfaces, output: fluxes
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)
235 integer :: n
236
237 do n = 1, nc
238 ! x-fluxes
239 call flux_koren_1d(cc(:, n), v(:, n, 1), nc, ngc)
240
241 ! y-fluxes (use temporary variables for efficiency)
242 cc_1d = cc(n, :)
243 v_1d = v(n, :, 2)
244 call flux_koren_1d(cc_1d, v_1d, nc, ngc)
245 v(n, :, 2) = v_1d ! Copy result
246 end do
247 end subroutine flux_koren_2d
248
249 !> Reconstruct "left" and "right" values at cell faces from cell-centered
250 !> data. The left value is for right-going waves, and the right value for
251 !> left-going waves.
252 subroutine reconstruct_lr_1d(nc, ngc, n_vars, cc, u_l, u_r, limiter)
253 integer, intent(in) :: nc !< Number of cells
254 integer, intent(in) :: ngc !< Number of ghost cells
255 integer, intent(in) :: n_vars !< Number of variables
256 real(dp), intent(in) :: cc(1-ngc:nc+ngc, n_vars) !< Cell-centered values
257 !> Reconstructed "left" values at every interface
258 real(dp), intent(inout) :: u_l(1:nc+1, n_vars)
259 !> Reconstructed "right" values at every interface
260 real(dp), intent(inout) :: u_r(1:nc+1, n_vars)
261 integer, intent(in) :: limiter !< Which limiter to use
262 real(dp) :: slopes(0:nc+1, n_vars)
263
264 associate(a=>cc(1:nc+2, :) - cc(0:nc+1, :), &
265 b=>cc(0:nc+1, :) - cc(-1:nc, :))
266 ! Compute limited slopes at the cell centers
267 slopes = af_limiter_apply(a, b, limiter)
268
269 ! Reconstruct values on the faces
270 u_l(1:nc+1, :) = cc(0:nc, :) + 0.5_dp * slopes(0:nc, :) ! left
271
272 if (af_limiter_symmetric(limiter)) then
273 u_r(1:nc+1, :) = cc(1:nc+1, :) - 0.5_dp * slopes(1:nc+1, :) ! right
274 else
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, :) ! right
277 end if
278 end associate
279 end subroutine reconstruct_lr_1d
280
281 !> Reconstruct upwind values at cell faces from cell-centered data.
282 subroutine reconstruct_upwind_1d(nc, ngc, n_vars, cc, u_f, limiter, direction_positive)
283 integer, intent(in) :: nc !< Number of cells
284 integer, intent(in) :: ngc !< Number of ghost cells
285 integer, intent(in) :: n_vars !< Number of variables
286 real(dp), intent(in) :: cc(1-ngc:nc+ngc, n_vars) !< Cell-centered values
287 !> Reconstructed values at every interface
288 real(dp), intent(inout) :: u_f(1:nc+1, n_vars)
289 integer, intent(in) :: limiter !< Which limiter to use
290 !> True means positive flow (to the "right"), false to the left
291 logical, intent(in) :: direction_positive(nc+1, n_vars)
292
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)
298 elsewhere
299 u_f = cc(1:nc+1, :) - 0.5_dp * &
300 af_limiter_apply(b(2:nc+2, :), a(2:nc+2, :), limiter)
301 end where
302 end associate
303 end subroutine reconstruct_upwind_1d
304
305 !> Compute flux for the KT scheme
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)
316
317 flux = 0.5_dp * (flux_l + flux_r - spread(wmax, 2, n_vars) * (u_r - u_l))
318 end subroutine flux_kurganovtadmor_1d
319
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 !< Time step
324 integer, intent(in) :: n_vars !< Number of cell-centered variables
325 integer, intent(in) :: i_cc(n_vars) !< All cell-centered indices
326 integer, intent(in) :: n_vars_flux !< Number of variables with fluxes
327 integer, intent(in) :: i_cc_flux(n_vars_flux) !< Cell-centered indices of flux variables
328 integer, intent(in) :: i_flux(n_vars_flux) !< Indices of fluxes
329 integer, intent(in) :: s_deriv !< State to compute derivatives from
330 integer, intent(in) :: n_prev !< Number of previous states
331 integer, intent(in) :: s_prev(n_prev) !< Previous states
332 real(dp), intent(in) :: w_prev(n_prev) !< Weights of previous states
333 integer, intent(in) :: s_out !< Output time state
334 !> Method to include source terms
335 procedure(subr_source) :: add_source_box
336 !> Number of time steps restrictions
337 integer, intent(in) :: n_dt
338 !> Maximal allowed time steps
339 real(dp), intent(out) :: dt_lim(n_dt)
340 !> If present, only update where the mask is true
341 procedure(subr_box_mask), optional :: get_mask
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)
347
348 nc = tree%n_cell
349 rfac = 0.0_dp ! Prevent warnings in 3D
350
351 dt_lim = 1e100_dp
352
353 !$omp parallel private(lvl, n, m, id, IJK, dt_dr, i_var, rfac, iv, mask, my_dt, tmp) &
354 !$omp &reduction(min:dt_lim)
355 my_dt = 1e100_dp
356
357 do lvl = 1, tree%highest_lvl
358 !$omp do
359 do n = 1, size(tree%lvls(lvl)%leaves)
360 id = tree%lvls(lvl)%leaves(n)
361 dt_dr = dt/tree%boxes(id)%dr
362
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)
366 else
367 mask = .true.
368 end if
369
370 ! Set output state to the weighted sum of previous states
371 do i_var = 1, n_vars
372 iv = i_cc(i_var)
373 tmp = 0.0_dp
374
375 do m = 1, n_prev
376 tmp = tmp + w_prev(m) * cc(dtimes(1:nc), iv+s_prev(m))
377 end do
378
379 cc(dtimes(1:nc), iv+s_out) = tmp
380 end do
381
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)
385
386#if NDIM == 1
387 do kji_do(1, nc)
388 if (mask(ijk)) then
389 cc(i, i_cc_flux+s_out) = cc(i, i_cc_flux+s_out) + &
390 dt_dr(1) * &
391 (fc(i, 1, i_flux) - fc(i+1, 1, i_flux))
392 end if
393 end do; close_do
394#elif NDIM == 2
395 if (tree%coord_t == af_cyl) then
396 call af_cyl_flux_factors(tree%boxes(id), rfac)
397 do kji_do(1, nc)
398 if (mask(ijk)) then
399 cc(i, j, i_cc_flux+s_out) = cc(i, j, i_cc_flux+s_out) + &
400 dt_dr(1) * (&
401 rfac(1, i) * fc(i, j, 1, i_flux) - &
402 rfac(2, i) * fc(i+1, j, 1, i_flux)) &
403 + dt_dr(2) * &
404 (fc(i, j, 2, i_flux) - fc(i, j+1, 2, i_flux))
405 end if
406 end do; close_do
407 else
408 do kji_do(1, nc)
409 if (mask(ijk)) then
410 cc(i, j, i_cc_flux+s_out) = cc(i, j, i_cc_flux+s_out) + &
411 dt_dr(1) * &
412 (fc(i, j, 1, i_flux) - fc(i+1, j, 1, i_flux)) &
413 + dt_dr(2) * &
414 (fc(i, j, 2, i_flux) - fc(i, j+1, 2, i_flux))
415 end if
416 end do; close_do
417 end if
418#elif NDIM == 3
419 do kji_do(1, nc)
420 if (mask(ijk)) then
421 cc(i, j, k, i_cc_flux+s_out) = cc(i, j, k, i_cc_flux+s_out) + &
422 dt_dr(1) * &
423 (fc(i, j, k, 1, i_flux) - fc(i+1, j, k, 1, i_flux)) + &
424 dt_dr(2) * &
425 (fc(i, j, k, 2, i_flux) - fc(i, j+1, k, 2, i_flux)) + &
426 dt_dr(3) * &
427 (fc(i, j, k, 3, i_flux) - fc(i, j, k+1, 3, i_flux))
428 end if
429 end do; close_do
430#endif
431 end associate
432 end do
433 !$omp end do
434 end do
435 !$omp end parallel
436 end subroutine flux_update_densities
437
438 !> Compute generic finite volume flux with a second order MUSCL scheme
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)
442 use m_af_restrict
443 use m_af_core
444 type(af_t), intent(inout) :: tree
445 integer, intent(in) :: n_vars !< Number of variables
446 integer, intent(in) :: i_cc(n_vars) !< Cell-centered variables
447 integer, intent(in) :: s_deriv !< State to compute derivatives from
448 integer, intent(in) :: i_flux(n_vars) !< Flux variables
449 !> Maximal time step, assuming a CFL number of 1.0
450 real(dp), intent(out) :: dt_lim
451 !> Compute the maximum wave speed
452 procedure(subr_max_wavespeed) :: max_wavespeed
453 !> Compute the flux from primitive variables
454 procedure(subr_flux) :: flux_from_primitives
455 !> Other flux contributions
456 procedure(subr_flux_modify) :: flux_modify
457 !> Potentially modify line densities
458 procedure(subr_line_modify) :: line_modify
459 !> Convert conservative variables to primitive ones
460 procedure(subr_prim_cons) :: to_primitive
461 !> Convert primitive variables to conservative ones
462 procedure(subr_prim_cons) :: to_conservative
463 !> Type of slope limiter to use for flux calculation
464 integer, intent(in) :: limiter
465 real(dp) :: my_dt
466
467 integer :: lvl, i
468
469 ! Ensure ghost cells near refinement boundaries can be properly filled
470 call af_restrict_ref_boundary(tree, i_cc+s_deriv)
471
472 dt_lim = 1e100_dp
473 my_dt = 1e100_dp
474
475 !$omp parallel private(lvl, i, my_dt) reduction(min:dt_lim)
476 do lvl = 1, tree%highest_lvl
477 !$omp do
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)
484 end do
485 !$omp end do
486 end do
487 !$omp end parallel
488
489 ! Compute coarse fluxes from the fine ones at refinement boundaries
490 call af_consistent_fluxes(tree, i_flux)
491
492 end subroutine flux_generic_tree
493
494 !> Compute generic finite volume flux with a second order MUSCL scheme
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)
498 use m_af_types
500 type(af_t), intent(inout) :: tree
501 integer, intent(in) :: id !< Id of box
502 integer, intent(in) :: nc !< Number of cells
503 integer, intent(in) :: n_vars !< Number of variables
504 integer, intent(in) :: i_cc(n_vars) !< Cell-centered variables
505 integer, intent(in) :: s_deriv !< State to compute derivatives from
506 integer, intent(in) :: i_flux(n_vars) !< Flux variables
507 !> Maximal time step, assuming a CFL number of 1.0
508 real(dp), intent(out) :: dt_lim
509 !> Compute the maximum wave speed
510 procedure(subr_max_wavespeed) :: max_wavespeed
511 !> Compute the flux from primitive variables on cell faces
512 procedure(subr_flux) :: flux_from_primitives
513 !> Modify flux for other contributions
514 procedure(subr_flux_modify) :: flux_modify
515 !> Potentially modify line densities
516 procedure(subr_line_modify) :: line_modify
517 !> Convert conservative variables to primitive ones
518 procedure(subr_prim_cons) :: to_primitive
519 !> Convert primitive variables to conservative ones
520 procedure(subr_prim_cons) :: to_conservative
521 !> Type of slope limiter to use for flux calculation
522 integer, intent(in) :: limiter
523
524
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)
533#if NDIM > 1
534 integer :: i
535#endif
536#if NDIM > 2
537 integer :: j
538#endif
539
540 inv_dr = 1/tree%boxes(id)%dr
541
542 ! Get two layers of ghost cell data
543 call af_gc2_box(tree, id, i_cc+s_deriv, cc)
544
545 ! This will contain the sum of CFL-related conditions. For example, one can
546 ! write dt * (vx/dx + vy/dy) < 1. This sum will then contain (vx/dx +
547 ! vy/dy).
548 cfl_sum = 0.0_dp
549
550 ! Jannis: Below, there are function calls in the inner part of a loop. When
551 ! I did some benchmarks, it was not significantly slower than using a buffer
552 ! and fewer functions calls.
553
554 do flux_dim = 1, ndim
555#if NDIM > 2
556 do j = 1, nc
557#endif
558#if NDIM > 1
559 do i = 1, nc
560#endif
561 ! Extract cell-centered values along a line
562 select case (flux_dim)
563#if NDIM == 1
564 case (1)
565 cc_line = cc(:, :)
566#elif NDIM == 2
567 case (1)
568 cc_line = cc(:, i, :)
569 case (2)
570 cc_line = cc(i, :, :)
571#elif NDIM == 3
572 case (1)
573 cc_line = cc(:, i, j, :)
574 case (2)
575 cc_line = cc(i, :, j, :)
576 case (3)
577 cc_line = cc(i, j, :, :)
578#endif
579 end select
580
581#if NDIM == 2
582 line_ix = [i]
583#elif NDIM == 3
584 line_ix = [i, j]
585#endif
586
587 ! Optionally modify data, e.g. to take into account a boundary
588 call line_modify(nc+4, n_vars, cc_line, flux_dim, &
589 tree%boxes(id), line_ix, s_deriv)
590
591 ! Flux computation is based on primitive variables (this can be a
592 ! dummy procedure)
593 call to_primitive(nc+4, n_vars, cc_line)
594
595 ! Reconstruct to cell faces, getting a left and right state at
596 ! every face
597 call reconstruct_lr_1d(nc, 2, n_vars, cc_line, u_l, u_r, limiter)
598
599 ! Determine the maximal wave speed for the left and right state
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)
602
603 ! Compute the fluxes for the left and right state
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)
608
609 ! Convert back to conservative form
610 call to_conservative(nc+1, n_vars, u_l)
611 call to_conservative(nc+1, n_vars, u_r)
612
613 ! Get maximum of left/right wave speed
614 w_l = max(w_l, w_r)
615
616 ! Get maximal wave speed for each cell center
617 cfl_sum_line = max(w_l(1:nc), w_l(2:nc+1)) * inv_dr(ndim)
618
619 ! Combine left and right fluxes to obtain a single flux
620 call flux_kurganovtadmor_1d(nc+1, n_vars, flux_l, flux_r, &
621 u_l, u_r, w_l, flux)
622
623 ! Potentially add other flux components
624 call flux_modify(nc+1, n_vars, flux_dim, flux, &
625 tree%boxes(id), line_ix, s_deriv)
626
627 ! Store the computed fluxes
628 select case (flux_dim)
629#if NDIM == 1
630 case (1)
631 tree%boxes(id)%fc(:, flux_dim, i_flux) = flux
632 cfl_sum = cfl_sum + cfl_sum_line
633#elif NDIM == 2
634 case (1)
635 tree%boxes(id)%fc(:, i, flux_dim, i_flux) = flux
636 cfl_sum(:, i) = cfl_sum(:, i) + cfl_sum_line
637 case (2)
638 tree%boxes(id)%fc(i, :, flux_dim, i_flux) = flux
639 cfl_sum(i, :) = cfl_sum(i, :) + cfl_sum_line
640#elif NDIM == 3
641 case (1)
642 tree%boxes(id)%fc(:, i, j, flux_dim, i_flux) = flux
643 cfl_sum(:, i, j) = cfl_sum(:, i, j) + cfl_sum_line
644 case (2)
645 tree%boxes(id)%fc(i, :, j, flux_dim, i_flux) = flux
646 cfl_sum(i, :, j) = cfl_sum(i, :, j) + cfl_sum_line
647 case (3)
648 tree%boxes(id)%fc(i, j, :, flux_dim, i_flux) = flux
649 cfl_sum(i, j, :) = cfl_sum(i, j, :) + cfl_sum_line
650#endif
651 end select
652#if NDIM > 1
653 end do
654#endif
655#if NDIM > 2
656 end do
657#endif
658 end do
659
660 ! Determine maximal time step
661 dt_lim = 1/max(maxval(cfl_sum), 1e-100_dp)
662
663 end subroutine flux_generic_box
664
665 !> Compute upwind fluxes
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)
668 use m_af_restrict
669 use m_af_core
670 type(af_t), intent(inout) :: tree
671 integer, intent(in) :: n_vars !< Number of variables
672 integer, intent(in) :: i_cc(n_vars) !< Cell-centered variables
673 integer, intent(in) :: s_deriv !< State to compute derivatives from
674 integer, intent(in) :: i_flux(n_vars) !< Flux variables
675 !> Number of time steps restrictions (first is CFL)
676 integer, intent(in) :: n_dt
677 !> Maximal allowed time steps, assuming a CFL number of 1.0
678 real(dp), intent(out) :: dt_lim(n_dt)
679 !> Method to compute fluxes
680 procedure(subr_flux_upwind) :: flux_upwind
681 !> Method to get direction of flux (positive or negative)
682 procedure(subr_flux_dir) :: flux_direction
683 !> Potentially modify line densities
684 procedure(subr_line_modify) :: line_modify
685 !> Type of slope limiter to use for flux calculation
686 integer, intent(in) :: limiter
687 integer :: lvl, i
688 real(dp) :: my_dt(n_dt)
689
690 ! Ensure ghost cells near refinement boundaries can be properly filled
691 call af_restrict_ref_boundary(tree, i_cc+s_deriv)
692
693 dt_lim = 1e100_dp
694 my_dt = 1e100_dp
695
696 !$omp parallel private(lvl, i, my_dt) reduction(min:dt_lim)
697 do lvl = 1, tree%highest_lvl
698 !$omp do
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)
704 end do
705 !$omp end do
706 end do
707 !$omp end parallel
708
709 ! Compute coarse fluxes from the fine ones at refinement boundaries
710 call af_consistent_fluxes(tree, i_flux)
711
712 end subroutine flux_upwind_tree
713
714 !> Compute generic finite volume flux
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)
717 use m_af_types
719 type(af_t), intent(inout) :: tree
720 integer, intent(in) :: id !< Id of box
721 integer, intent(in) :: nc !< Number of cells
722 integer, intent(in) :: n_vars !< Number of variables
723 integer, intent(in) :: i_cc(n_vars) !< Cell-centered variables
724 integer, intent(in) :: s_deriv !< State to compute derivatives from
725 integer, intent(in) :: i_flux(n_vars) !< Flux variables
726 !> Number of time steps restrictions (first is CFL)
727 integer, intent(in) :: n_dt
728 real(dp), intent(inout) :: dt_lim(n_dt) !< Time step restrictions
729 !> Method to compute fluxes
730 procedure(subr_flux_upwind) :: flux_upwind
731 !> Method to get direction of flux (positive or negative)
732 procedure(subr_flux_dir) :: flux_direction
733 !> Potentially modify line densities
734 procedure(subr_line_modify) :: line_modify
735 !> Type of slope limiter to use for flux calculation
736 integer, intent(in) :: limiter
737
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)
746#if NDIM > 1
747 integer :: i
748#endif
749#if NDIM > 2
750 integer :: j
751#endif
752
753 ! Get two layers of ghost cell data
754 call af_gc2_box(tree, id, i_cc+s_deriv, cc)
755
756 ! This will contain the sum of CFL-related conditions. For example, one can
757 ! write dt * (vx/dx + vy/dy) < 1. This sum will then contain (vx/dx +
758 ! vy/dy).
759 cfl_sum = 0.0_dp
760 dt_lim(2:) = 1e100_dp
761 other_dt = 1e100_dp
762
763 associate(fc => tree%boxes(id)%fc)
764 do flux_dim = 1, ndim
765#if NDIM > 2
766 do j = 1, nc
767#endif
768#if NDIM > 1
769 do i = 1, nc
770#endif
771 ! Extract cell-centered values along a line
772 select case (flux_dim)
773#if NDIM == 1
774 case (1)
775 cc_line = cc(:, :)
776#elif NDIM == 2
777 case (1)
778 cc_line = cc(:, i, :)
779 case (2)
780 cc_line = cc(i, :, :)
781#elif NDIM == 3
782 case (1)
783 cc_line = cc(:, i, j, :)
784 case (2)
785 cc_line = cc(i, :, j, :)
786 case (3)
787 cc_line = cc(i, j, :, :)
788#endif
789 end select
790
791#if NDIM == 2
792 line_ix = [i]
793#elif NDIM == 3
794 line_ix = [i, j]
795#endif
796 call flux_direction(tree%boxes(id), line_ix, s_deriv, &
797 n_vars, flux_dim, direction_positive)
798
799 ! Optionally modify data, e.g. to take into account a boundary
800 call line_modify(nc+4, n_vars, cc_line, flux_dim, &
801 tree%boxes(id), line_ix, s_deriv)
802
803 ! Reconstruct to cell faces
804 call reconstruct_upwind_1d(nc, 2, n_vars, cc_line, u_l, &
805 limiter, direction_positive)
806
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)
810
811 ! Store the computed fluxes
812 select case (flux_dim)
813#if NDIM == 1
814 case (1)
815 fc(:, flux_dim, i_flux) = flux
816 cfl_sum = cfl_sum + cfl_sum_line
817#elif NDIM == 2
818 case (1)
819 fc(:, i, flux_dim, i_flux) = flux
820 cfl_sum(:, i) = cfl_sum(:, i) + cfl_sum_line
821 case (2)
822 fc(i, :, flux_dim, i_flux) = flux
823 cfl_sum(i, :) = cfl_sum(i, :) + cfl_sum_line
824#elif NDIM == 3
825 case (1)
826 fc(:, i, j, flux_dim, i_flux) = flux
827 cfl_sum(:, i, j) = cfl_sum(:, i, j) + cfl_sum_line
828 case (2)
829 fc(i, :, j, flux_dim, i_flux) = flux
830 cfl_sum(i, :, j) = cfl_sum(i, :, j) + cfl_sum_line
831 case (3)
832 fc(i, j, :, flux_dim, i_flux) = flux
833 cfl_sum(i, j, :) = cfl_sum(i, j, :) + cfl_sum_line
834#endif
835 end select
836#if NDIM > 1
837 end do
838#endif
839#if NDIM > 2
840 end do
841#endif
842 end do
843 end associate
844
845 ! Determine maximal CFL time step
846 dt_lim(1) = 1/maxval(cfl_sum)
847
848 end subroutine flux_upwind_box
849
850 !> Extract cell-centered data along a line in a box, including a single layer
851 !> of ghost cells. This is convenient to get extra variables in a flux
852 !> computation.
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(:) !< Indices of the variables
856 integer, intent(in) :: flux_dim !< Dimension of flux computation
857 integer, intent(in) :: line_ix(ndim-1) !< Index of line
858 real(dp), intent(inout) :: cc_line(box%n_cell+2, size(ivs))
859
860 select case (flux_dim)
861#if NDIM == 1
862 case (1)
863 cc_line = box%cc(:, ivs)
864#elif NDIM == 2
865 case (1)
866 cc_line = box%cc(:, line_ix(1), ivs)
867 case (2)
868 cc_line = box%cc(line_ix(1), :, ivs)
869#elif NDIM == 3
870 case (1)
871 cc_line = box%cc(:, line_ix(1), line_ix(2), ivs)
872 case (2)
873 cc_line = box%cc(line_ix(1), :, line_ix(2), ivs)
874 case (3)
875 cc_line = box%cc(line_ix(1), line_ix(2), :, ivs)
876#endif
877 end select
878 end subroutine flux_get_line_cc
879
880 !> Extract cell-centered data along a line in a box, including a single layer
881 !> of ghost cells. This is convenient to get extra variables in a flux
882 !> computation.
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 !< Index of the variable
886 integer, intent(in) :: flux_dim !< Dimension of flux computation
887 integer, intent(in) :: line_ix(ndim-1) !< Index of line
888 real(dp), intent(inout) :: cc_line(box%n_cell+2)
889
890 select case (flux_dim)
891#if NDIM == 1
892 case (1)
893 cc_line = box%cc(:, iv)
894#elif NDIM == 2
895 case (1)
896 cc_line = box%cc(:, line_ix(1), iv)
897 case (2)
898 cc_line = box%cc(line_ix(1), :, iv)
899#elif NDIM == 3
900 case (1)
901 cc_line = box%cc(:, line_ix(1), line_ix(2), iv)
902 case (2)
903 cc_line = box%cc(line_ix(1), :, line_ix(2), iv)
904 case (3)
905 cc_line = box%cc(line_ix(1), line_ix(2), :, iv)
906#endif
907 end select
908 end subroutine flux_get_line_1cc
909
910 !> Extract face-centered data along a line in a box. This is convenient to get
911 !> extra variables in a flux computation.
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(:) !< Indices of the variables
915 integer, intent(in) :: flux_dim !< Dimension of flux computation
916 integer, intent(in) :: line_ix(ndim-1) !< Index of line
917 real(dp), intent(inout) :: fc_line(box%n_cell+1, size(ivs))
918
919 select case (flux_dim)
920#if NDIM == 1
921 case (1)
922 fc_line = box%fc(:, flux_dim, ivs)
923#elif NDIM == 2
924 case (1)
925 fc_line = box%fc(:, line_ix(1), flux_dim, ivs)
926 case (2)
927 fc_line = box%fc(line_ix(1), :, flux_dim, ivs)
928#elif NDIM == 3
929 case (1)
930 fc_line = box%fc(:, line_ix(1), line_ix(2), flux_dim, ivs)
931 case (2)
932 fc_line = box%fc(line_ix(1), :, line_ix(2), flux_dim, ivs)
933 case (3)
934 fc_line = box%fc(line_ix(1), line_ix(2), :, flux_dim, ivs)
935#endif
936 end select
937 end subroutine flux_get_line_fc
938
939 !> Extract face-centered data along a line in a box. This is convenient to get
940 !> extra variables in a flux computation.
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 !< Index of the variable
944 integer, intent(in) :: flux_dim !< Dimension of flux computation
945 integer, intent(in) :: line_ix(ndim-1) !< Index of line
946 real(dp), intent(inout) :: fc_line(box%n_cell+1)
947
948 select case (flux_dim)
949#if NDIM == 1
950 case (1)
951 fc_line = box%fc(:, flux_dim, iv)
952#elif NDIM == 2
953 case (1)
954 fc_line = box%fc(:, line_ix(1), flux_dim, iv)
955 case (2)
956 fc_line = box%fc(line_ix(1), :, flux_dim, iv)
957#elif NDIM == 3
958 case (1)
959 fc_line = box%fc(:, line_ix(1), line_ix(2), flux_dim, iv)
960 case (2)
961 fc_line = box%fc(line_ix(1), :, line_ix(2), flux_dim, iv)
962 case (3)
963 fc_line = box%fc(line_ix(1), line_ix(2), :, flux_dim, iv)
964#endif
965 end select
966 end subroutine flux_get_line_1fc
967
968 !> Compute flux according to Koren limiter
969 subroutine flux_koren_3d(cc, v, nc, ngc)
970 !> Number of cells
971 integer, intent(in) :: nc
972 !> Number of ghost cells
973 integer, intent(in) :: ngc
974 !> Cell-centered values
975 real(dp), intent(in) :: cc(1-ngc:nc+ngc, 1-ngc:nc+ngc, 1-ngc:nc+ngc)
976 !> Input: velocities at interfaces, output: fluxes
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)
979 integer :: n, m
980
981 do m = 1, nc
982 do n = 1, nc
983 ! x-fluxes
984 call flux_koren_1d(cc(:, n, m), v(:, n, m, 1), &
985 nc, ngc)
986
987 ! y-fluxes (use temporary variables for efficiency)
988 cc_1d = cc(n, :, m)
989 v_1d = v(n, :, m, 2)
990 call flux_koren_1d(cc_1d, v_1d, nc, ngc)
991 v(n, :, m, 2) = v_1d ! Copy result
992
993 ! z-fluxes (use temporary variables for efficiency)
994 cc_1d = cc(n, m, :)
995 v_1d = v(n, m, :, 3)
996 call flux_koren_1d(cc_1d, v_1d, nc, ngc)
997 v(n, m, :, 3) = v_1d ! Copy result
998 end do
999 end do
1000 end subroutine flux_koren_3d
1001
1002 !> Compute flux with first order upwind scheme
1003 subroutine flux_upwind_1d(cc, v, nc, ngc)
1004 integer, intent(in) :: nc !< Number of cells
1005 integer, intent(in) :: ngc !< Number of ghost cells
1006 real(dp), intent(in) :: cc(1-ngc:nc+ngc) !< Cell-centered values
1007 !> Input: velocities at interfaces, output: fluxes
1008 real(dp), intent(inout) :: v(1:nc+1)
1009 integer :: n
1010
1011 do n = 1, nc+1
1012 if (v(n) < 0.0_dp) then
1013 v(n) = v(n) * cc(n)
1014 else ! v(n) > 0
1015 v(n) = v(n) * cc(n-1)
1016 end if
1017 end do
1018
1019 end subroutine flux_upwind_1d
1020
1021 !> Compute flux with first order upwind scheme
1022 subroutine flux_upwind_2d(cc, v, nc, ngc)
1023 integer, intent(in) :: nc !< Number of cells
1024 integer, intent(in) :: ngc !< Number of ghost cells
1025 !> Cell-centered values
1026 real(dp), intent(in) :: cc(1-ngc:nc+ngc, 1-ngc:nc+ngc)
1027 !> Input: velocities at interfaces, output: fluxes
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)
1030 integer :: n
1031
1032 do n = 1, nc
1033 ! x-fluxes
1034 call flux_upwind_1d(cc(:, n), v(:, n, 1), nc, ngc)
1035
1036 ! y-fluxes (use temporary variables for efficiency)
1037 cc_1d = cc(n, :)
1038 v_1d = v(n, :, 2)
1039 call flux_upwind_1d(cc_1d, v_1d, nc, ngc)
1040 v(n, :, 2) = v_1d ! Copy result
1041 end do
1042 end subroutine flux_upwind_2d
1043
1044 !> Compute flux with first order upwind scheme
1045 subroutine flux_upwind_3d(cc, v, nc, ngc)
1046 !> Number of cells
1047 integer, intent(in) :: nc
1048 !> Number of ghost cells
1049 integer, intent(in) :: ngc
1050 !> Cell-centered values
1051 real(dp), intent(in) :: cc(1-ngc:nc+ngc, 1-ngc:nc+ngc, 1-ngc:nc+ngc)
1052 !> Input: velocities at interfaces, output: fluxes
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)
1055 integer :: n, m
1056
1057 do m = 1, nc
1058 do n = 1, nc
1059 ! x-fluxes
1060 call flux_upwind_1d(cc(:, n, m), v(:, n, m, 1), &
1061 nc, ngc)
1062
1063 ! y-fluxes (use temporary variables for efficiency)
1064 cc_1d = cc(n, :, m)
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 ! Copy result
1068
1069 ! z-fluxes (use temporary variables for efficiency)
1070 cc_1d = cc(n, m, :)
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 ! Copy result
1074 end do
1075 end do
1076 end subroutine flux_upwind_3d
1077
1078 !> Dummy conversion between primitive and conservative variables
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
1083
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
1095
1096 subroutine flux_dummy_modify(nf, n_var, flux_dim, flux, box, line_ix, s_deriv)
1097 integer, intent(in) :: nf !< Number of cell faces
1098 integer, intent(in) :: n_var !< Number of variables
1099 integer, intent(in) :: flux_dim !< In which dimension fluxes are computed
1100 real(dp), intent(inout) :: flux(nf, n_var) !< Flux that will be modified
1101 type(box_t), intent(in) :: box !< Current box
1102 integer, intent(in) :: line_ix(ndim-1) !< Index of line for dim /= flux_dim
1103 integer, intent(in) :: s_deriv !< State to compute derivatives from
1104 end subroutine flux_dummy_modify
1105
1106 subroutine flux_dummy_line_modify(n_cc, n_var, cc_line, flux_dim, box, &
1107 line_ix, s_deriv)
1108 integer, intent(in) :: n_cc !< Number of cell centers
1109 integer, intent(in) :: n_var !< Number of variables
1110 real(dp), intent(inout) :: cc_line(n_cc, n_var) !< Flux that will be modified
1111 integer, intent(in) :: flux_dim !< In which dimension fluxes are computed
1112 type(box_t), intent(in) :: box !< Current box
1113 integer, intent(in) :: line_ix(ndim-1) !< Index of line for dim /= flux_dim
1114 integer, intent(in) :: s_deriv !< State to compute derivatives from
1115 end subroutine flux_dummy_line_modify
1116
1117end module m_af_flux_schemes
This module contains the core routines of Afivo, namely those that deal with initializing and changin...
Definition: m_af_core.f90:3
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...
Definition: m_af_types.f90:3
Type which stores all the boxes and levels, as well as some information about the number of boxes,...
Definition: m_af_types.f90:326
The basic building block of afivo: a box with cell-centered and face centered data,...
Definition: m_af_types.f90:286