Afivo
0.3
|
Module containing a couple flux schemes for solving hyperbolic problems explicitly, as well as handling diffusion explicitly. More...
Data Types | |
interface | subr_prim_cons |
interface | subr_max_wavespeed |
interface | subr_flux |
interface | subr_flux_upwind |
interface | subr_flux_modify |
interface | subr_flux_dir |
interface | subr_line_modify |
interface | subr_source |
interface | subr_box_mask |
Functions/Subroutines | |
subroutine, public | flux_diff_1d (cc, dc, inv_dr, nc, ngc) |
Compute diffusive flux (second order) More... | |
subroutine, public | flux_diff_2d (cc, dc, inv_dr, nc, ngc) |
Compute diffusive flux (second order) More... | |
subroutine, public | flux_diff_3d (cc, dc, inv_dr, nc, ngc) |
Compute diffusive flux (second order) More... | |
subroutine, public | flux_koren_1d (cc, v, nc, ngc) |
Compute flux according to Koren limiter. More... | |
subroutine, public | flux_koren_2d (cc, v, nc, ngc) |
Compute flux according to Koren limiter. More... | |
subroutine, public | reconstruct_lr_1d (nc, ngc, n_vars, cc, u_l, u_r, limiter) |
Reconstruct "left" and "right" values at cell faces from cell-centered data. The left value is for right-going waves, and the right value for left-going waves. More... | |
subroutine | reconstruct_upwind_1d (nc, ngc, n_vars, cc, u_f, limiter, direction_positive) |
Reconstruct upwind values at cell faces from cell-centered data. More... | |
subroutine, public | flux_kurganovtadmor_1d (n_values, n_vars, flux_l, flux_r, u_l, u_r, wmax, flux) |
Compute flux for the KT scheme. More... | |
subroutine, public | flux_update_densities (tree, dt, n_vars, i_cc, n_vars_flux, i_cc_flux, i_flux, s_deriv, n_prev, s_prev, w_prev, s_out, add_source_box, n_dt, dt_lim, get_mask) |
subroutine, public | flux_generic_tree (tree, n_vars, i_cc, s_deriv, i_flux, dt_lim, max_wavespeed, flux_from_primitives, flux_modify, line_modify, to_primitive, to_conservative, limiter) |
Compute generic finite volume flux with a second order MUSCL scheme. More... | |
subroutine, public | flux_generic_box (tree, id, nc, n_vars, i_cc, s_deriv, i_flux, dt_lim, max_wavespeed, flux_from_primitives, flux_modify, line_modify, to_primitive, to_conservative, limiter) |
Compute generic finite volume flux with a second order MUSCL scheme. More... | |
subroutine, public | flux_upwind_tree (tree, n_vars, i_cc, s_deriv, i_flux, n_dt, dt_lim, flux_upwind, flux_direction, line_modify, limiter) |
Compute upwind fluxes. More... | |
subroutine, public | flux_upwind_box (tree, id, nc, n_vars, i_cc, s_deriv, i_flux, n_dt, dt_lim, flux_upwind, flux_direction, line_modify, limiter) |
Compute generic finite volume flux. More... | |
subroutine, public | flux_get_line_cc (box, ivs, flux_dim, line_ix, cc_line) |
Extract cell-centered data along a line in a box, including a single layer of ghost cells. This is convenient to get extra variables in a flux computation. More... | |
subroutine, public | flux_get_line_1cc (box, iv, flux_dim, line_ix, cc_line) |
Extract cell-centered data along a line in a box, including a single layer of ghost cells. This is convenient to get extra variables in a flux computation. More... | |
subroutine, public | flux_get_line_fc (box, ivs, flux_dim, line_ix, fc_line) |
Extract face-centered data along a line in a box. This is convenient to get extra variables in a flux computation. More... | |
subroutine, public | flux_get_line_1fc (box, iv, flux_dim, line_ix, fc_line) |
Extract face-centered data along a line in a box. This is convenient to get extra variables in a flux computation. More... | |
subroutine, public | flux_koren_3d (cc, v, nc, ngc) |
Compute flux according to Koren limiter. More... | |
subroutine, public | flux_upwind_1d (cc, v, nc, ngc) |
Compute flux with first order upwind scheme. More... | |
subroutine, public | flux_upwind_2d (cc, v, nc, ngc) |
Compute flux with first order upwind scheme. More... | |
subroutine, public | flux_upwind_3d (cc, v, nc, ngc) |
Compute flux with first order upwind scheme. More... | |
subroutine, public | flux_dummy_conversion (n_values, n_vars, u) |
Dummy conversion between primitive and conservative variables. More... | |
subroutine, public | flux_dummy_source (box, dt, n_vars, i_cc, s_deriv, s_out, n_dt, dt_lim, mask) |
subroutine, public | flux_dummy_modify (nf, n_var, flux_dim, flux, box, line_ix, s_deriv) |
subroutine, public | flux_dummy_line_modify (n_cc, n_var, cc_line, flux_dim, box, line_ix, s_deriv) |
Module containing a couple flux schemes for solving hyperbolic problems explicitly, as well as handling diffusion explicitly.
subroutine, public m_af_flux_schemes::flux_diff_1d | ( | real(dp), dimension(1-ngc:nc+ngc), intent(in) | cc, |
real(dp), dimension(1:nc+1), intent(inout) | dc, | ||
real(dp), intent(in) | inv_dr, | ||
integer, intent(in) | nc, | ||
integer, intent(in) | ngc | ||
) |
Compute diffusive flux (second order)
[in] | nc | Number of cells |
[in] | ngc | Number of ghost cells |
[in] | cc | Cell-centered values |
[in,out] | dc | Input: diffusion coeff. at interfaces, output: fluxes |
[in] | inv_dr | Inverse grid spacing |
Definition at line 129 of file m_af_flux_schemes.f90.
subroutine, public m_af_flux_schemes::flux_diff_2d | ( | real(dp), dimension(1-ngc:nc+ngc, 1-ngc:nc+ngc), intent(in) | cc, |
real(dp), dimension(1:nc+1, 1:nc+1, 2), intent(inout) | dc, | ||
real(dp), dimension(2), intent(in) | inv_dr, | ||
integer, intent(in) | nc, | ||
integer, intent(in) | ngc | ||
) |
Compute diffusive flux (second order)
[in] | nc | Number of cells |
[in] | ngc | Number of ghost cells |
[in] | cc | Cell-centered values |
[in,out] | dc | Input: diffusion coeff. at interfaces, output: fluxes |
[in] | inv_dr | Inverse grid spacing |
Definition at line 144 of file m_af_flux_schemes.f90.
subroutine, public m_af_flux_schemes::flux_diff_3d | ( | real(dp), dimension(1-ngc:nc+ngc, 1-ngc:nc+ngc, 1-ngc:nc+ngc), intent(in) | cc, |
real(dp), dimension(1:nc+1, 1:nc+1, 1:nc+1, 3), intent(inout) | dc, | ||
real(dp), dimension(3), intent(in) | inv_dr, | ||
integer, intent(in) | nc, | ||
integer, intent(in) | ngc | ||
) |
Compute diffusive flux (second order)
[in] | nc | Number of cells |
[in] | ngc | Number of ghost cells |
[in] | cc | Cell-centered values |
[in,out] | dc | Input: diffusion coeff. at interfaces, output: fluxes |
[in] | inv_dr | Inverse grid spacing |
Definition at line 168 of file m_af_flux_schemes.f90.
subroutine, public m_af_flux_schemes::flux_koren_1d | ( | real(dp), dimension(1-ngc:nc+ngc), intent(in) | cc, |
real(dp), dimension(1:nc+1), intent(inout) | v, | ||
integer, intent(in) | nc, | ||
integer, intent(in) | ngc | ||
) |
Compute flux according to Koren limiter.
[in] | nc | Number of cells |
[in] | ngc | Number of ghost cells |
[in] | cc | Cell-centered values |
[in,out] | v | Input: velocities at interfaces, output: fluxes |
Definition at line 204 of file m_af_flux_schemes.f90.
subroutine, public m_af_flux_schemes::flux_koren_2d | ( | real(dp), dimension(1-ngc:nc+ngc, 1-ngc:nc+ngc), intent(in) | cc, |
real(dp), dimension(1:nc+1, 1:nc+1, 2), intent(inout) | v, | ||
integer, intent(in) | nc, | ||
integer, intent(in) | ngc | ||
) |
Compute flux according to Koren limiter.
[in] | nc | Number of cells |
[in] | ngc | Number of ghost cells |
[in] | cc | Cell-centered values |
[in,out] | v | Input: velocities at interfaces, output: fluxes |
Definition at line 227 of file m_af_flux_schemes.f90.
subroutine, public m_af_flux_schemes::reconstruct_lr_1d | ( | integer, intent(in) | nc, |
integer, intent(in) | ngc, | ||
integer, intent(in) | n_vars, | ||
real(dp), dimension(1-ngc:nc+ngc, n_vars), intent(in) | cc, | ||
real(dp), dimension(1:nc+1, n_vars), intent(inout) | u_l, | ||
real(dp), dimension(1:nc+1, n_vars), intent(inout) | u_r, | ||
integer, intent(in) | limiter | ||
) |
Reconstruct "left" and "right" values at cell faces from cell-centered data. The left value is for right-going waves, and the right value for left-going waves.
[in] | nc | Number of cells |
[in] | ngc | Number of ghost cells |
[in] | n_vars | Number of variables |
[in] | cc | Cell-centered values |
[in,out] | u_l | Reconstructed "left" values at every interface |
[in,out] | u_r | Reconstructed "right" values at every interface |
[in] | limiter | Which limiter to use |
Definition at line 252 of file m_af_flux_schemes.f90.
|
private |
Reconstruct upwind values at cell faces from cell-centered data.
[in] | nc | Number of cells |
[in] | ngc | Number of ghost cells |
[in] | n_vars | Number of variables |
[in] | cc | Cell-centered values |
[in,out] | u_f | Reconstructed values at every interface |
[in] | limiter | Which limiter to use |
[in] | direction_positive | True means positive flow (to the "right"), false to the left |
Definition at line 282 of file m_af_flux_schemes.f90.
subroutine, public m_af_flux_schemes::flux_kurganovtadmor_1d | ( | integer, intent(in) | n_values, |
integer, intent(in) | n_vars, | ||
real(dp), dimension(n_values, n_vars), intent(in) | flux_l, | ||
real(dp), dimension(n_values, n_vars), intent(in) | flux_r, | ||
real(dp), dimension(n_values, n_vars), intent(in) | u_l, | ||
real(dp), dimension(n_values, n_vars), intent(in) | u_r, | ||
real(dp), dimension(n_values), intent(in) | wmax, | ||
real(dp), dimension(n_values, n_vars), intent(inout) | flux | ||
) |
Compute flux for the KT scheme.
Definition at line 306 of file m_af_flux_schemes.f90.
subroutine, public m_af_flux_schemes::flux_update_densities | ( | type(af_t), intent(inout) | tree, |
real(dp), intent(in) | dt, | ||
integer, intent(in) | n_vars, | ||
integer, dimension(n_vars), intent(in) | i_cc, | ||
integer, intent(in) | n_vars_flux, | ||
integer, dimension(n_vars_flux), intent(in) | i_cc_flux, | ||
integer, dimension(n_vars_flux), intent(in) | i_flux, | ||
integer, intent(in) | s_deriv, | ||
integer, intent(in) | n_prev, | ||
integer, dimension(n_prev), intent(in) | s_prev, | ||
real(dp), dimension(n_prev), intent(in) | w_prev, | ||
integer, intent(in) | s_out, | ||
procedure(subr_source) | add_source_box, | ||
integer, intent(in) | n_dt, | ||
real(dp), dimension(n_dt), intent(out) | dt_lim, | ||
procedure(subr_box_mask), optional | get_mask | ||
) |
[in] | dt | Time step |
[in] | n_vars | Number of cell-centered variables |
[in] | i_cc | All cell-centered indices |
[in] | n_vars_flux | Number of variables with fluxes |
[in] | i_cc_flux | Cell-centered indices of flux variables |
[in] | i_flux | Indices of fluxes |
[in] | s_deriv | State to compute derivatives from |
[in] | n_prev | Number of previous states |
[in] | s_prev | Previous states |
[in] | w_prev | Weights of previous states |
[in] | s_out | Output time state |
add_source_box | Method to include source terms | |
[in] | n_dt | Number of time steps restrictions |
[out] | dt_lim | Maximal allowed time steps |
get_mask | If present, only update where the mask is true |
Definition at line 320 of file m_af_flux_schemes.f90.
subroutine, public m_af_flux_schemes::flux_generic_tree | ( | type(af_t), intent(inout) | tree, |
integer, intent(in) | n_vars, | ||
integer, dimension(n_vars), intent(in) | i_cc, | ||
integer, intent(in) | s_deriv, | ||
integer, dimension(n_vars), intent(in) | i_flux, | ||
real(dp), intent(out) | dt_lim, | ||
procedure(subr_max_wavespeed) | max_wavespeed, | ||
procedure(subr_flux) | flux_from_primitives, | ||
procedure(subr_flux_modify) | flux_modify, | ||
procedure(subr_line_modify) | line_modify, | ||
procedure(subr_prim_cons) | to_primitive, | ||
procedure(subr_prim_cons) | to_conservative, | ||
integer, intent(in) | limiter | ||
) |
Compute generic finite volume flux with a second order MUSCL scheme.
[in] | n_vars | Number of variables |
[in] | i_cc | Cell-centered variables |
[in] | s_deriv | State to compute derivatives from |
[in] | i_flux | Flux variables |
[out] | dt_lim | Maximal time step, assuming a CFL number of 1.0 |
max_wavespeed | Compute the maximum wave speed | |
flux_from_primitives | Compute the flux from primitive variables | |
flux_modify | Other flux contributions | |
line_modify | Potentially modify line densities | |
to_primitive | Convert conservative variables to primitive ones | |
to_conservative | Convert primitive variables to conservative ones | |
[in] | limiter | Type of slope limiter to use for flux calculation |
Definition at line 439 of file m_af_flux_schemes.f90.
subroutine, public m_af_flux_schemes::flux_generic_box | ( | type(af_t), intent(inout) | tree, |
integer, intent(in) | id, | ||
integer, intent(in) | nc, | ||
integer, intent(in) | n_vars, | ||
integer, dimension(n_vars), intent(in) | i_cc, | ||
integer, intent(in) | s_deriv, | ||
integer, dimension(n_vars), intent(in) | i_flux, | ||
real(dp), intent(out) | dt_lim, | ||
procedure(subr_max_wavespeed) | max_wavespeed, | ||
procedure(subr_flux) | flux_from_primitives, | ||
procedure(subr_flux_modify) | flux_modify, | ||
procedure(subr_line_modify) | line_modify, | ||
procedure(subr_prim_cons) | to_primitive, | ||
procedure(subr_prim_cons) | to_conservative, | ||
integer, intent(in) | limiter | ||
) |
Compute generic finite volume flux with a second order MUSCL scheme.
[in] | id | Id of box |
[in] | nc | Number of cells |
[in] | n_vars | Number of variables |
[in] | i_cc | Cell-centered variables |
[in] | s_deriv | State to compute derivatives from |
[in] | i_flux | Flux variables |
[out] | dt_lim | Maximal time step, assuming a CFL number of 1.0 |
max_wavespeed | Compute the maximum wave speed | |
flux_from_primitives | Compute the flux from primitive variables on cell faces | |
flux_modify | Modify flux for other contributions | |
line_modify | Potentially modify line densities | |
to_primitive | Convert conservative variables to primitive ones | |
to_conservative | Convert primitive variables to conservative ones | |
[in] | limiter | Type of slope limiter to use for flux calculation |
Definition at line 495 of file m_af_flux_schemes.f90.
subroutine, public m_af_flux_schemes::flux_upwind_tree | ( | type(af_t), intent(inout) | tree, |
integer, intent(in) | n_vars, | ||
integer, dimension(n_vars), intent(in) | i_cc, | ||
integer, intent(in) | s_deriv, | ||
integer, dimension(n_vars), intent(in) | i_flux, | ||
integer, intent(in) | n_dt, | ||
real(dp), dimension(n_dt), intent(out) | dt_lim, | ||
procedure(subr_flux_upwind) | flux_upwind, | ||
procedure(subr_flux_dir) | flux_direction, | ||
procedure(subr_line_modify) | line_modify, | ||
integer, intent(in) | limiter | ||
) |
Compute upwind fluxes.
[in] | n_vars | Number of variables |
[in] | i_cc | Cell-centered variables |
[in] | s_deriv | State to compute derivatives from |
[in] | i_flux | Flux variables |
[in] | n_dt | Number of time steps restrictions (first is CFL) |
[out] | dt_lim | Maximal allowed time steps, assuming a CFL number of 1.0 |
flux_upwind | Method to compute fluxes | |
flux_direction | Method to get direction of flux (positive or negative) | |
line_modify | Potentially modify line densities | |
[in] | limiter | Type of slope limiter to use for flux calculation |
Definition at line 666 of file m_af_flux_schemes.f90.
subroutine, public m_af_flux_schemes::flux_upwind_box | ( | type(af_t), intent(inout) | tree, |
integer, intent(in) | id, | ||
integer, intent(in) | nc, | ||
integer, intent(in) | n_vars, | ||
integer, dimension(n_vars), intent(in) | i_cc, | ||
integer, intent(in) | s_deriv, | ||
integer, dimension(n_vars), intent(in) | i_flux, | ||
integer, intent(in) | n_dt, | ||
real(dp), dimension(n_dt), intent(inout) | dt_lim, | ||
procedure(subr_flux_upwind) | flux_upwind, | ||
procedure(subr_flux_dir) | flux_direction, | ||
procedure(subr_line_modify) | line_modify, | ||
integer, intent(in) | limiter | ||
) |
Compute generic finite volume flux.
[in] | id | Id of box |
[in] | nc | Number of cells |
[in] | n_vars | Number of variables |
[in] | i_cc | Cell-centered variables |
[in] | s_deriv | State to compute derivatives from |
[in] | i_flux | Flux variables |
[in] | n_dt | Number of time steps restrictions (first is CFL) |
[in,out] | dt_lim | Time step restrictions |
flux_upwind | Method to compute fluxes | |
flux_direction | Method to get direction of flux (positive or negative) | |
line_modify | Potentially modify line densities | |
[in] | limiter | Type of slope limiter to use for flux calculation |
Definition at line 715 of file m_af_flux_schemes.f90.
subroutine, public m_af_flux_schemes::flux_get_line_cc | ( | type(box_t), intent(in) | box, |
integer, dimension(:), intent(in) | ivs, | ||
integer, intent(in) | flux_dim, | ||
integer, dimension(ndim-1), intent(in) | line_ix, | ||
real(dp), dimension(box%n_cell+2, size(ivs)), intent(inout) | cc_line | ||
) |
Extract cell-centered data along a line in a box, including a single layer of ghost cells. This is convenient to get extra variables in a flux computation.
[in] | ivs | Indices of the variables |
[in] | flux_dim | Dimension of flux computation |
[in] | line_ix | Index of line |
Definition at line 853 of file m_af_flux_schemes.f90.
subroutine, public m_af_flux_schemes::flux_get_line_1cc | ( | type(box_t), intent(in) | box, |
integer, intent(in) | iv, | ||
integer, intent(in) | flux_dim, | ||
integer, dimension(ndim-1), intent(in) | line_ix, | ||
real(dp), dimension(box%n_cell+2), intent(inout) | cc_line | ||
) |
Extract cell-centered data along a line in a box, including a single layer of ghost cells. This is convenient to get extra variables in a flux computation.
[in] | iv | Index of the variable |
[in] | flux_dim | Dimension of flux computation |
[in] | line_ix | Index of line |
Definition at line 883 of file m_af_flux_schemes.f90.
subroutine, public m_af_flux_schemes::flux_get_line_fc | ( | type(box_t), intent(in) | box, |
integer, dimension(:), intent(in) | ivs, | ||
integer, intent(in) | flux_dim, | ||
integer, dimension(ndim-1), intent(in) | line_ix, | ||
real(dp), dimension(box%n_cell+1, size(ivs)), intent(inout) | fc_line | ||
) |
Extract face-centered data along a line in a box. This is convenient to get extra variables in a flux computation.
[in] | ivs | Indices of the variables |
[in] | flux_dim | Dimension of flux computation |
[in] | line_ix | Index of line |
Definition at line 912 of file m_af_flux_schemes.f90.
subroutine, public m_af_flux_schemes::flux_get_line_1fc | ( | type(box_t), intent(in) | box, |
integer, intent(in) | iv, | ||
integer, intent(in) | flux_dim, | ||
integer, dimension(ndim-1), intent(in) | line_ix, | ||
real(dp), dimension(box%n_cell+1), intent(inout) | fc_line | ||
) |
Extract face-centered data along a line in a box. This is convenient to get extra variables in a flux computation.
[in] | iv | Index of the variable |
[in] | flux_dim | Dimension of flux computation |
[in] | line_ix | Index of line |
Definition at line 941 of file m_af_flux_schemes.f90.
subroutine, public m_af_flux_schemes::flux_koren_3d | ( | real(dp), dimension(1-ngc:nc+ngc, 1-ngc:nc+ngc, 1-ngc:nc+ngc), intent(in) | cc, |
real(dp), dimension(1:nc+1, 1:nc+1, 1:nc+1, 3), intent(inout) | v, | ||
integer, intent(in) | nc, | ||
integer, intent(in) | ngc | ||
) |
Compute flux according to Koren limiter.
[in] | nc | Number of cells |
[in] | ngc | Number of ghost cells |
[in] | cc | Cell-centered values |
[in,out] | v | Input: velocities at interfaces, output: fluxes |
Definition at line 969 of file m_af_flux_schemes.f90.
subroutine, public m_af_flux_schemes::flux_upwind_1d | ( | real(dp), dimension(1-ngc:nc+ngc), intent(in) | cc, |
real(dp), dimension(1:nc+1), intent(inout) | v, | ||
integer, intent(in) | nc, | ||
integer, intent(in) | ngc | ||
) |
Compute flux with first order upwind scheme.
[in] | nc | Number of cells |
[in] | ngc | Number of ghost cells |
[in] | cc | Cell-centered values |
[in,out] | v | Input: velocities at interfaces, output: fluxes |
Definition at line 1003 of file m_af_flux_schemes.f90.
subroutine, public m_af_flux_schemes::flux_upwind_2d | ( | real(dp), dimension(1-ngc:nc+ngc, 1-ngc:nc+ngc), intent(in) | cc, |
real(dp), dimension(1:nc+1, 1:nc+1, 2), intent(inout) | v, | ||
integer, intent(in) | nc, | ||
integer, intent(in) | ngc | ||
) |
Compute flux with first order upwind scheme.
[in] | nc | Number of cells |
[in] | ngc | Number of ghost cells |
[in] | cc | Cell-centered values |
[in,out] | v | Input: velocities at interfaces, output: fluxes |
Definition at line 1022 of file m_af_flux_schemes.f90.
subroutine, public m_af_flux_schemes::flux_upwind_3d | ( | real(dp), dimension(1-ngc:nc+ngc, 1-ngc:nc+ngc, 1-ngc:nc+ngc), intent(in) | cc, |
real(dp), dimension(1:nc+1, 1:nc+1, 1:nc+1, 3), intent(inout) | v, | ||
integer, intent(in) | nc, | ||
integer, intent(in) | ngc | ||
) |
Compute flux with first order upwind scheme.
[in] | nc | Number of cells |
[in] | ngc | Number of ghost cells |
[in] | cc | Cell-centered values |
[in,out] | v | Input: velocities at interfaces, output: fluxes |
Definition at line 1045 of file m_af_flux_schemes.f90.
subroutine, public m_af_flux_schemes::flux_dummy_conversion | ( | integer, intent(in) | n_values, |
integer, intent(in) | n_vars, | ||
real(dp), dimension(n_values, n_vars), intent(inout) | u | ||
) |
Dummy conversion between primitive and conservative variables.
Definition at line 1079 of file m_af_flux_schemes.f90.
subroutine, public m_af_flux_schemes::flux_dummy_source | ( | type(box_t), intent(inout) | box, |
real(dp), intent(in) | dt, | ||
integer, intent(in) | n_vars, | ||
integer, dimension(n_vars), intent(in) | i_cc, | ||
integer, intent(in) | s_deriv, | ||
integer, intent(in) | s_out, | ||
integer, intent(in) | n_dt, | ||
real(dp), dimension(n_dt), intent(inout) | dt_lim, | ||
logical, dimension(dtimes(box%n_cell)), intent(in) | mask | ||
) |
Definition at line 1084 of file m_af_flux_schemes.f90.
subroutine, public m_af_flux_schemes::flux_dummy_modify | ( | integer, intent(in) | nf, |
integer, intent(in) | n_var, | ||
integer, intent(in) | flux_dim, | ||
real(dp), dimension(nf, n_var), intent(inout) | flux, | ||
type(box_t), intent(in) | box, | ||
integer, dimension(ndim-1), intent(in) | line_ix, | ||
integer, intent(in) | s_deriv | ||
) |
[in] | nf | Number of cell faces |
[in] | n_var | Number of variables |
[in] | flux_dim | In which dimension fluxes are computed |
[in,out] | flux | Flux that will be modified |
[in] | box | Current box |
[in] | line_ix | Index of line for dim /= flux_dim |
[in] | s_deriv | State to compute derivatives from |
Definition at line 1096 of file m_af_flux_schemes.f90.
subroutine, public m_af_flux_schemes::flux_dummy_line_modify | ( | integer, intent(in) | n_cc, |
integer, intent(in) | n_var, | ||
real(dp), dimension(n_cc, n_var), intent(inout) | cc_line, | ||
integer, intent(in) | flux_dim, | ||
type(box_t), intent(in) | box, | ||
integer, dimension(ndim-1), intent(in) | line_ix, | ||
integer, intent(in) | s_deriv | ||
) |
[in] | n_cc | Number of cell centers |
[in] | n_var | Number of variables |
[in,out] | cc_line | Flux that will be modified |
[in] | flux_dim | In which dimension fluxes are computed |
[in] | box | Current box |
[in] | line_ix | Index of line for dim /= flux_dim |
[in] | s_deriv | State to compute derivatives from |
Definition at line 1106 of file m_af_flux_schemes.f90.