Afivo  0.3
Data Types | Functions/Subroutines
m_af_flux_schemes Module Reference

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)
 

Detailed Description

Module containing a couple flux schemes for solving hyperbolic problems explicitly, as well as handling diffusion explicitly.

Function/Subroutine Documentation

◆ flux_diff_1d()

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)

Parameters
[in]ncNumber of cells
[in]ngcNumber of ghost cells
[in]ccCell-centered values
[in,out]dcInput: diffusion coeff. at interfaces, output: fluxes
[in]inv_drInverse grid spacing

Definition at line 129 of file m_af_flux_schemes.f90.

◆ flux_diff_2d()

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)

Parameters
[in]ncNumber of cells
[in]ngcNumber of ghost cells
[in]ccCell-centered values
[in,out]dcInput: diffusion coeff. at interfaces, output: fluxes
[in]inv_drInverse grid spacing

Definition at line 144 of file m_af_flux_schemes.f90.

◆ flux_diff_3d()

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)

Parameters
[in]ncNumber of cells
[in]ngcNumber of ghost cells
[in]ccCell-centered values
[in,out]dcInput: diffusion coeff. at interfaces, output: fluxes
[in]inv_drInverse grid spacing

Definition at line 168 of file m_af_flux_schemes.f90.

◆ flux_koren_1d()

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.

Parameters
[in]ncNumber of cells
[in]ngcNumber of ghost cells
[in]ccCell-centered values
[in,out]vInput: velocities at interfaces, output: fluxes

Definition at line 204 of file m_af_flux_schemes.f90.

◆ flux_koren_2d()

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.

Parameters
[in]ncNumber of cells
[in]ngcNumber of ghost cells
[in]ccCell-centered values
[in,out]vInput: velocities at interfaces, output: fluxes

Definition at line 227 of file m_af_flux_schemes.f90.

◆ reconstruct_lr_1d()

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.

Parameters
[in]ncNumber of cells
[in]ngcNumber of ghost cells
[in]n_varsNumber of variables
[in]ccCell-centered values
[in,out]u_lReconstructed "left" values at every interface
[in,out]u_rReconstructed "right" values at every interface
[in]limiterWhich limiter to use

Definition at line 252 of file m_af_flux_schemes.f90.

◆ reconstruct_upwind_1d()

subroutine m_af_flux_schemes::reconstruct_upwind_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_f,
integer, intent(in)  limiter,
logical, dimension(nc+1, n_vars), intent(in)  direction_positive 
)
private

Reconstruct upwind values at cell faces from cell-centered data.

Parameters
[in]ncNumber of cells
[in]ngcNumber of ghost cells
[in]n_varsNumber of variables
[in]ccCell-centered values
[in,out]u_fReconstructed values at every interface
[in]limiterWhich limiter to use
[in]direction_positiveTrue means positive flow (to the "right"), false to the left

Definition at line 282 of file m_af_flux_schemes.f90.

◆ flux_kurganovtadmor_1d()

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.

◆ flux_update_densities()

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 
)
Parameters
[in]dtTime step
[in]n_varsNumber of cell-centered variables
[in]i_ccAll cell-centered indices
[in]n_vars_fluxNumber of variables with fluxes
[in]i_cc_fluxCell-centered indices of flux variables
[in]i_fluxIndices of fluxes
[in]s_derivState to compute derivatives from
[in]n_prevNumber of previous states
[in]s_prevPrevious states
[in]w_prevWeights of previous states
[in]s_outOutput time state
add_source_boxMethod to include source terms
[in]n_dtNumber of time steps restrictions
[out]dt_limMaximal allowed time steps
get_maskIf present, only update where the mask is true

Definition at line 320 of file m_af_flux_schemes.f90.

◆ flux_generic_tree()

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.

Parameters
[in]n_varsNumber of variables
[in]i_ccCell-centered variables
[in]s_derivState to compute derivatives from
[in]i_fluxFlux variables
[out]dt_limMaximal time step, assuming a CFL number of 1.0
max_wavespeedCompute the maximum wave speed
flux_from_primitivesCompute the flux from primitive variables
flux_modifyOther flux contributions
line_modifyPotentially modify line densities
to_primitiveConvert conservative variables to primitive ones
to_conservativeConvert primitive variables to conservative ones
[in]limiterType of slope limiter to use for flux calculation

Definition at line 434 of file m_af_flux_schemes.f90.

◆ flux_generic_box()

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.

Parameters
[in]idId of box
[in]ncNumber of cells
[in]n_varsNumber of variables
[in]i_ccCell-centered variables
[in]s_derivState to compute derivatives from
[in]i_fluxFlux variables
[out]dt_limMaximal time step, assuming a CFL number of 1.0
max_wavespeedCompute the maximum wave speed
flux_from_primitivesCompute the flux from primitive variables on cell faces
flux_modifyModify flux for other contributions
line_modifyPotentially modify line densities
to_primitiveConvert conservative variables to primitive ones
to_conservativeConvert primitive variables to conservative ones
[in]limiterType of slope limiter to use for flux calculation

Definition at line 490 of file m_af_flux_schemes.f90.

◆ flux_upwind_tree()

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.

Parameters
[in]n_varsNumber of variables
[in]i_ccCell-centered variables
[in]s_derivState to compute derivatives from
[in]i_fluxFlux variables
[in]n_dtNumber of time steps restrictions (first is CFL)
[out]dt_limMaximal allowed time steps, assuming a CFL number of 1.0
flux_upwindMethod to compute fluxes
flux_directionMethod to get direction of flux (positive or negative)
line_modifyPotentially modify line densities
[in]limiterType of slope limiter to use for flux calculation

Definition at line 661 of file m_af_flux_schemes.f90.

◆ flux_upwind_box()

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.

Parameters
[in]idId of box
[in]ncNumber of cells
[in]n_varsNumber of variables
[in]i_ccCell-centered variables
[in]s_derivState to compute derivatives from
[in]i_fluxFlux variables
[in]n_dtNumber of time steps restrictions (first is CFL)
[in,out]dt_limTime step restrictions
flux_upwindMethod to compute fluxes
flux_directionMethod to get direction of flux (positive or negative)
line_modifyPotentially modify line densities
[in]limiterType of slope limiter to use for flux calculation

Definition at line 710 of file m_af_flux_schemes.f90.

◆ flux_get_line_cc()

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.

Parameters
[in]ivsIndices of the variables
[in]flux_dimDimension of flux computation
[in]line_ixIndex of line

Definition at line 848 of file m_af_flux_schemes.f90.

◆ flux_get_line_1cc()

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.

Parameters
[in]ivIndex of the variable
[in]flux_dimDimension of flux computation
[in]line_ixIndex of line

Definition at line 878 of file m_af_flux_schemes.f90.

◆ flux_get_line_fc()

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.

Parameters
[in]ivsIndices of the variables
[in]flux_dimDimension of flux computation
[in]line_ixIndex of line

Definition at line 907 of file m_af_flux_schemes.f90.

◆ flux_get_line_1fc()

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.

Parameters
[in]ivIndex of the variable
[in]flux_dimDimension of flux computation
[in]line_ixIndex of line

Definition at line 936 of file m_af_flux_schemes.f90.

◆ flux_koren_3d()

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.

Parameters
[in]ncNumber of cells
[in]ngcNumber of ghost cells
[in]ccCell-centered values
[in,out]vInput: velocities at interfaces, output: fluxes

Definition at line 964 of file m_af_flux_schemes.f90.

◆ flux_upwind_1d()

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.

Parameters
[in]ncNumber of cells
[in]ngcNumber of ghost cells
[in]ccCell-centered values
[in,out]vInput: velocities at interfaces, output: fluxes

Definition at line 998 of file m_af_flux_schemes.f90.

◆ flux_upwind_2d()

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.

Parameters
[in]ncNumber of cells
[in]ngcNumber of ghost cells
[in]ccCell-centered values
[in,out]vInput: velocities at interfaces, output: fluxes

Definition at line 1017 of file m_af_flux_schemes.f90.

◆ flux_upwind_3d()

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.

Parameters
[in]ncNumber of cells
[in]ngcNumber of ghost cells
[in]ccCell-centered values
[in,out]vInput: velocities at interfaces, output: fluxes

Definition at line 1040 of file m_af_flux_schemes.f90.

◆ flux_dummy_conversion()

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 1074 of file m_af_flux_schemes.f90.

◆ flux_dummy_source()

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 1079 of file m_af_flux_schemes.f90.

◆ flux_dummy_modify()

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 
)
Parameters
[in]nfNumber of cell faces
[in]n_varNumber of variables
[in]flux_dimIn which dimension fluxes are computed
[in,out]fluxFlux that will be modified
[in]boxCurrent box
[in]line_ixIndex of line for dim /= flux_dim
[in]s_derivState to compute derivatives from

Definition at line 1091 of file m_af_flux_schemes.f90.

◆ flux_dummy_line_modify()

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 
)
Parameters
[in]n_ccNumber of cell centers
[in]n_varNumber of variables
[in,out]cc_lineFlux that will be modified
[in]flux_dimIn which dimension fluxes are computed
[in]boxCurrent box
[in]line_ixIndex of line for dim /= flux_dim
[in]s_derivState to compute derivatives from

Definition at line 1101 of file m_af_flux_schemes.f90.