IMEX_SfloW2D  0.9
Shallowwatergranularflowmodel
Functions/Subroutines | Variables
solver_2d Module Reference

Numerical solver. More...

Functions/Subroutines

subroutine allocate_solver_variables
 Memory allocation. More...
 
subroutine deallocate_solver_variables
 Memory deallocation. More...
 
subroutine check_solve
 Masking of cells to solve. More...
 
subroutine timestep
 Time-step computation. More...
 
subroutine timestep2
 Time-step computation. More...
 
subroutine imex_rk_solver
 Runge-Kutta integration. More...
 
subroutine solve_rk_step (Bj, Bprimej_x, Bprimej_y, grav3_surf, qj, qj_old, a_tilde, a_dirk, a_diag)
 Runge-Kutta single step integration. More...
 
subroutine lnsrch (Bj, Bprimej_x, Bprimej_y, grav3_surf, qj_rel_NR_old, qj_org, qj_old, scal_f_old, grad_f, desc_dir, coeff_f, qj_rel, scal_f, right_term, stpmax, check)
 Search the descent stepsize. More...
 
subroutine eval_f (Bj, Bprimej_x, Bprimej_y, grav3_surf, qj, qj_old, a_tilde, a_dirk, a_diag, coeff_f, f_nl, scal_f)
 Evaluate the nonlinear system. More...
 
subroutine eval_jacobian (Bj, Bprimej_x, Bprimej_y, grav3_surf, qj_rel, qj_org, coeff_f, left_matrix)
 Evaluate the jacobian. More...
 
subroutine eval_explicit_terms (q_expl, expl_terms)
 Evaluate the explicit terms. More...
 
subroutine eval_hyperbolic_terms (q_expl, divFlux)
 Semidiscrete finite volume central scheme. More...
 
subroutine eval_flux_kt
 Semidiscrete numerical fluxes. More...
 
subroutine average_kt (a1, a2, w1, w2, w_avg)
 averaged KT flux More...
 
subroutine eval_flux_gforce
 Numerical fluxes GFORCE. More...
 
subroutine eval_flux_lxf
 Numerical fluxes Lax-Friedrichs. More...
 
subroutine reconstruction
 Linear reconstruction. More...
 
subroutine eval_speeds
 Characteristic speeds. More...
 
subroutine limit (v, z, limiter, slope_lim)
 Slope limiter. More...
 
real *8 function minmod (a, b)
 
real *8 function maxmod (a, b)
 

Variables

real *8, dimension(:,:,:), allocatable q
 Conservative variables. More...
 
real *8, dimension(:,:,:), allocatable q0
 Conservative variables. More...
 
real *8, dimension(:,:,:), allocatable q_fv
 Solution of the finite-volume semidiscrete cheme. More...
 
real *8, dimension(:,:,:), allocatable q_interfacel
 Reconstructed value at the left of the x-interface. More...
 
real *8, dimension(:,:,:), allocatable q_interfacer
 Reconstructed value at the right of the x-interface. More...
 
real *8, dimension(:,:,:), allocatable q_interfaceb
 Reconstructed value at the bottom of the y-interface. More...
 
real *8, dimension(:,:,:), allocatable q_interfacet
 Reconstructed value at the top of the y-interface. More...
 
real *8, dimension(:,:,:), allocatable a_interface_xneg
 Local speeds at the left of the x-interface. More...
 
real *8, dimension(:,:,:), allocatable a_interface_xpos
 Local speeds at the right of the x-interface. More...
 
real *8, dimension(:,:,:), allocatable a_interface_yneg
 Local speeds at the bottom of the y-interface. More...
 
real *8, dimension(:,:,:), allocatable a_interface_ypos
 Local speeds at the top of the y-interface. More...
 
real *8, dimension(:,:,:), allocatable h_interface_x
 Semidiscrete numerical interface fluxes. More...
 
real *8, dimension(:,:,:), allocatable h_interface_y
 Semidiscrete numerical interface fluxes. More...
 
real *8, dimension(:,:,:), allocatable qp
 Physical variables ( $\alpha_1, p_1, p_2, \rho u, w, T$) More...
 
real *8, dimension(:,:), allocatable source_xy
 Array defining fraction of cells affected by source term. More...
 
logical, dimension(:,:), allocatable solve_mask
 
logical, dimension(:,:), allocatable solve_mask0
 
real *8 dt
 Time step. More...
 
logical, dimension(:,:), allocatable mask22
 
logical, dimension(:,:), allocatable mask21
 
logical, dimension(:,:), allocatable mask11
 
logical, dimension(:,:), allocatable mask12
 
integer i_rk
 loop counter for the RK iteration More...
 
real *8, dimension(:,:), allocatable a_tilde_ij
 Butcher Tableau for the explicit part of the Runge-Kutta scheme. More...
 
real *8, dimension(:,:), allocatable a_dirk_ij
 Butcher Tableau for the implicit part of the Runge-Kutta scheme. More...
 
real *8, dimension(:), allocatable omega_tilde
 Coefficients for the explicit part of the Runge-Kutta scheme. More...
 
real *8, dimension(:), allocatable omega
 Coefficients for the implicit part of the Runge-Kutta scheme. More...
 
real *8, dimension(:), allocatable a_tilde
 Explicit coeff. for the hyperbolic part for a single step of the R-K scheme. More...
 
real *8, dimension(:), allocatable a_dirk
 Explicit coeff. for the non-hyp. part for a single step of the R-K scheme. More...
 
real *8 a_diag
 Implicit coeff. for the non-hyp. part for a single step of the R-K scheme. More...
 
real *8, dimension(:,:,:,:), allocatable q_rk
 Intermediate solutions of the Runge-Kutta scheme. More...
 
real *8, dimension(:,:,:,:), allocatable divflux
 Intermediate hyperbolic terms of the Runge-Kutta scheme. More...
 
real *8, dimension(:,:,:,:), allocatable nh
 Intermediate non-hyperbolic terms of the Runge-Kutta scheme. More...
 
real *8, dimension(:,:,:,:), allocatable si_nh
 Intermediate semi-implicit non-hyperbolic terms of the Runge-Kutta scheme. More...
 
real *8, dimension(:,:,:,:), allocatable expl_terms
 Intermediate explicit terms of the Runge-Kutta scheme. More...
 
real *8, dimension(:,:), allocatable divfluxj
 Local Intermediate hyperbolic terms of the Runge-Kutta scheme. More...
 
real *8, dimension(:,:), allocatable nhj
 Local Intermediate non-hyperbolic terms of the Runge-Kutta scheme. More...
 
real *8, dimension(:,:), allocatable expl_terms_j
 Local Intermediate explicit terms of the Runge-Kutta scheme. More...
 
real *8, dimension(:,:), allocatable si_nhj
 Local Intermediate semi-impl non-hyperbolic terms of the Runge-Kutta scheme. More...
 
logical normalize_q
 Flag for the normalization of the array q in the implicit solution scheme. More...
 
logical normalize_f
 Flag for the normalization of the array f in the implicit solution scheme. More...
 
logical opt_search_nl
 Flag for the search of optimal step size in the implicit solution scheme. More...
 
real *8, dimension(:,:,:), allocatable residual_term
 Sum of all the terms of the equations except the transient term. More...
 
real *8 t_imex1
 
real *8 t_imex2
 

Detailed Description

Numerical solver.

This module contains the variables and the subroutines for the numerical solution of the equations.

Date
07/10/2016
Author
Mattia de' Michieli Vitturi

Function/Subroutine Documentation

subroutine solver_2d::allocate_solver_variables ( )

Memory allocation.

This subroutine allocate the memory for the variables of the solver module.

Date
07/10/2016
Author
Mattia de' Michieli Vitturi

Definition at line 141 of file solver_2d.f90.

Here is the caller graph for this function:

subroutine solver_2d::average_kt ( real*8, dimension(:), intent(in)  a1,
real*8, dimension(:), intent(in)  a2,
real*8, dimension(:), intent(in)  w1,
real*8, dimension(:), intent(in)  w2,
real*8, dimension(:), intent(out)  w_avg 
)

averaged KT flux

This subroutine compute n averaged flux from the fluxes at the two sides of a cell interface and the max an min speed at the two sides.

Parameters
[in]aLspeed at one side of the interface
[in]aRspeed at the other side of the interface
[in]wLfluxes at one side of the interface
[in]wRfluxes at the other side of the interface
[out]w_avgarray of averaged fluxes
Date
07/10/2016
Author
Mattia de' Michieli Vitturi

Definition at line 2007 of file solver_2d.f90.

Here is the caller graph for this function:

subroutine solver_2d::check_solve ( )

Masking of cells to solve.

This subroutine compute a 2D array of logicals defining the cells where the systems of equations have to be solved. It is defined according to the positive thickness in the cell and in the neighbour cells

Date
20/04/2017
Author
Mattia de' Michieli Vitturi

Definition at line 404 of file solver_2d.f90.

subroutine solver_2d::deallocate_solver_variables ( )

Memory deallocation.

This subroutine de-allocate the memory for the variables of the solver module.

Date
07/10/2016
Author
Mattia de' Michieli Vitturi

Definition at line 338 of file solver_2d.f90.

Here is the caller graph for this function:

subroutine solver_2d::eval_explicit_terms ( real*8, dimension(n_vars,comp_cells_x,comp_cells_y), intent(in)  q_expl,
real*8, dimension(n_eqns,comp_cells_x,comp_cells_y), intent(out)  expl_terms 
)

Evaluate the explicit terms.

This subroutine evaluate the explicit terms (non-fluxes) of the non-linear system with respect to the conservative variables.

Parameters
[in]q_explconservative variables
[out]expl_termsexplicit terms
Date
07/10/2016
Author
Mattia de' Michieli Vitturi

Definition at line 1677 of file solver_2d.f90.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine solver_2d::eval_f ( real*8, intent(in)  Bj,
real*8, intent(in)  Bprimej_x,
real*8, intent(in)  Bprimej_y,
real*8, intent(in)  grav3_surf,
real*8, dimension(n_vars), intent(in)  qj,
real*8, dimension(n_vars), intent(in)  qj_old,
real*8, dimension(n_rk), intent(in)  a_tilde,
real*8, dimension(n_rk), intent(in)  a_dirk,
real*8, intent(in)  a_diag,
real*8, dimension(n_eqns), intent(in)  coeff_f,
real*8, dimension(n_eqns), intent(out)  f_nl,
real*8, intent(out)  scal_f 
)

Evaluate the nonlinear system.

This subroutine evaluate the value of the nonlinear system in the state defined by the variables qj.

Parameters
[in]Bjtopography at the cell center
[in]Bprimejtopography slope at the cell center
[in]qjconservative variables
[in]qj_oldconservative variables at the old time step
[in]a_tildeexplicit coefficients for the hyperbolic terms
[in]a_dirkexplicit coefficients for the non-hyperbolic terms
[in]a_diagimplicit coefficient for the non-hyperbolic term
[in]coeff_fcoefficient to rescale the nonlinear functions
[out]f_nlvalues of the nonlinear functions
[out]scal_fvalue of the scalar function f=0.5*<F,F>
Date
07/10/2016
Author
Mattia de' Michieli Vitturi

Definition at line 1537 of file solver_2d.f90.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine solver_2d::eval_flux_gforce ( )

Numerical fluxes GFORCE.

Date
07/10/2016
Author
Mattia de' Michieli Vitturi

Definition at line 2044 of file solver_2d.f90.

Here is the caller graph for this function:

subroutine solver_2d::eval_flux_kt ( )

Semidiscrete numerical fluxes.

This subroutine evaluates the numerical fluxes H at the cells interfaces according to Kurganov et al. 2001.

Author
Mattia de' Michieli Vitturi
Date
16/08/2011

Definition at line 1864 of file solver_2d.f90.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine solver_2d::eval_flux_lxf ( )

Numerical fluxes Lax-Friedrichs.

Date
07/10/2016
Author
Mattia de' Michieli Vitturi

Definition at line 2058 of file solver_2d.f90.

Here is the caller graph for this function:

subroutine solver_2d::eval_hyperbolic_terms ( real*8, dimension(n_vars,comp_cells_x,comp_cells_y), intent(in)  q_expl,
real*8, dimension(n_eqns,comp_cells_x,comp_cells_y), intent(out)  divFlux 
)

Semidiscrete finite volume central scheme.

This subroutine compute the divergence part of the system of the eqns, with a modified version of the finite volume scheme from Kurganov et al. 2001, where the reconstruction at the cells interfaces is applied to a set of physical variables derived from the conservative vriables.

Parameters
[in]q_explconservative variables
[out]divFluxdivergence term
Date
07/10/2016
Author
Mattia de' Michieli Vitturi

Definition at line 1724 of file solver_2d.f90.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine solver_2d::eval_jacobian ( real*8, intent(in)  Bj,
real*8, intent(in)  Bprimej_x,
real*8, intent(in)  Bprimej_y,
real*8, intent(in)  grav3_surf,
real*8, dimension(n_vars), intent(in)  qj_rel,
real*8, dimension(n_vars), intent(in)  qj_org,
real*8, dimension(n_eqns), intent(in)  coeff_f,
real*8, dimension(n_eqns,n_vars), intent(out)  left_matrix 
)

Evaluate the jacobian.

This subroutine evaluate the jacobian of the non-linear system with respect to the conservative variables.

Parameters
[in]Bjtopography at the cell center
[in]Bprimej_xtopography x-slope at the cell center
[in]Bprimej_ytopography y-slope at the cell center
[in]grav3_surf
[in]qj_relrelative variation (qj=qj_rel*qj_org)
[in]qj_orgconservative variables at the old time step
[in]coeff_fcoefficient to rescale the nonlinear functions
[out]left_matrixmatrix from the linearization of the system
Date
07/10/2016
Author
Mattia de' Michieli Vitturi

Definition at line 1604 of file solver_2d.f90.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine solver_2d::eval_speeds ( )

Characteristic speeds.

This subroutine evaluates the largest characteristic speed at the cells interfaces from the reconstructed states.

Author
Mattia de' Michieli Vitturi
Date
16/08/2011

Definition at line 2523 of file solver_2d.f90.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine solver_2d::imex_rk_solver ( )

Runge-Kutta integration.

This subroutine integrate the hyperbolic conservation law with non-hyperbolic terms using an implicit-explicit runge-kutta scheme. The fluxes are integrated explicitely while the non-hyperbolic terms are integrated implicitely.

Date
07/10/2016
Author
Mattia de' Michieli Vitturi

Definition at line 598 of file solver_2d.f90.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine solver_2d::limit ( real*8, dimension(3), intent(in)  v,
real*8, dimension(3), intent(in)  z,
integer, intent(in)  limiter,
real*8, intent(out)  slope_lim 
)

Slope limiter.

This subroutine limits the slope of the linear reconstruction of the physical variables, accordingly to the parameter "solve_limiter":

  • 'none' => no limiter (constant value);
  • 'minmod' => minmod slope;
  • 'superbee' => superbee limiter (Roe, 1985);
  • 'van_leer' => monotonized central-difference limiter (van Leer, 1977)
Parameters
[in]v3-point stencil value array
[in]z3-point stencil location array
[in]limiterinteger defining the limiter choice
[out]slope_limlimited slope
Date
07/10/2016
Author
Mattia de' Michieli Vitturi

Definition at line 2609 of file solver_2d.f90.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine solver_2d::lnsrch ( real*8, intent(in)  Bj,
real*8, intent(in)  Bprimej_x,
real*8, intent(in)  Bprimej_y,
real*8, intent(in)  grav3_surf,
real*8, dimension(:), intent(in)  qj_rel_NR_old,
real*8, dimension(:), intent(in)  qj_org,
real*8, dimension(:), intent(in)  qj_old,
real*8, intent(in)  scal_f_old,
real*8, dimension(:), intent(in)  grad_f,
real*8, dimension(:), intent(inout)  desc_dir,
real*8, dimension(:), intent(in)  coeff_f,
real*8, dimension(:), intent(out)  qj_rel,
real*8, intent(out)  scal_f,
real*8, dimension(n_eqns), intent(out)  right_term,
real*8, intent(in)  stpmax,
logical, intent(out)  check 
)

Search the descent stepsize.

This subroutine search for the lenght of the descent step in order to have a decrease in the nonlinear function.

Parameters
[in]Bjtopography at the cell center
[in]Bprimej_xtopography x-slope at the cell center
[in]Bprimej_ytopography y-slope at the cell center
[in]qj_rel_NR_old
[in]qj_org
[in]qj_old
[in]scal_f_old
[in]grad_f
[in,out]desc_dir
[in]coeff_f
[out]qj_rel
[out]scal_f
[out]right_term
[in]stpmax
[out]check
Date
07/10/2016
Author
Mattia de' Michieli Vitturi
Parameters
[in]qj_rel_nr_oldInitial point
[in]qj_orgInitial point
[in]qj_oldInitial point
[in]grad_fGradient at xold
[in]scal_f_oldValue of the function at xold
[in,out]desc_dirDescent direction (usually Newton direction)
[in]coeff_fCoefficients to rescale the nonlinear function
[out]qj_relUpdated solution
[out]scal_fValue of the scalar function at x
[out]right_termValue of the scalar function at x
[out]checkOutput quantity check is false on a normal exit

Definition at line 1313 of file solver_2d.f90.

Here is the call graph for this function:

Here is the caller graph for this function:

real*8 function solver_2d::maxmod ( real*8  a,
real*8  b 
)

Definition at line 2687 of file solver_2d.f90.

Here is the caller graph for this function:

real*8 function solver_2d::minmod ( real*8  a,
real*8  b 
)

Definition at line 2666 of file solver_2d.f90.

Here is the caller graph for this function:

subroutine solver_2d::reconstruction ( )

Linear reconstruction.

In this subroutine a linear reconstruction with slope limiters is applied to a set of variables describing the state of the system, according to the input parameter reconstr_variables.

Author
Mattia de' Michieli Vitturi
Date
15/08/2011

Definition at line 2077 of file solver_2d.f90.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine solver_2d::solve_rk_step ( real*8, intent(in)  Bj,
real*8, intent(in)  Bprimej_x,
real*8, intent(in)  Bprimej_y,
real*8, intent(in)  grav3_surf,
real*8, dimension(n_vars), intent(inout)  qj,
real*8, dimension(n_vars), intent(in)  qj_old,
real*8, dimension(n_rk), intent(in)  a_tilde,
real*8, dimension(n_rk), intent(in)  a_dirk,
real*8, intent(in)  a_diag 
)

Runge-Kutta single step integration.

This subroutine find the solution of the non-linear system given the a step of the implicit-explicit Runge-Kutta scheme for a cell:
$ Q^{(i)} = Q^n - dt \sum_{j=1}^{i-1}\tilde{a}_{j}\partial_x F(Q^{(j)}) + dt \sum_{j=1}^{i-1} a_j NH(Q^{(j)}) + dt a_{diag} NH(Q^{(i)}) $

Parameters
[in]Bjtopography at the cell center
[in]Bprimejtopography slope at the cell center
[in,out]qjconservative variables
[in]qj_oldconservative variables at the old time step
[in]a_tildeexplicit coefficents for the fluxes
[in]a_dirkexplicit coefficient for the non-hyperbolic terms
[in]a_diagimplicit coefficient for the non-hyperbolic terms
Date
07/10/2016
Author
Mattia de' Michieli Vitturi

Definition at line 987 of file solver_2d.f90.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine solver_2d::timestep ( )

Time-step computation.

This subroutine evaluate the maximum time step according to the CFL condition. The local speed are evaluated with the characteristic polynomial of the Jacobian of the fluxes.

Date
07/10/2016
Author
Mattia de' Michieli Vitturi

Definition at line 448 of file solver_2d.f90.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine solver_2d::timestep2 ( )

Time-step computation.

This subroutine evaluate the maximum time step according to the CFL condition. The local speed are evaluated with the characteristic polynomial of the Jacobian of the fluxes.

Date
07/10/2016
Author
Mattia de' Michieli Vitturi

Definition at line 525 of file solver_2d.f90.

Here is the call graph for this function:

Here is the caller graph for this function:

Variable Documentation

real*8 solver_2d::a_diag

Implicit coeff. for the non-hyp. part for a single step of the R-K scheme.

Definition at line 88 of file solver_2d.f90.

real*8, dimension(:), allocatable solver_2d::a_dirk

Explicit coeff. for the non-hyp. part for a single step of the R-K scheme.

Definition at line 86 of file solver_2d.f90.

real*8, dimension(:,:), allocatable solver_2d::a_dirk_ij

Butcher Tableau for the implicit part of the Runge-Kutta scheme.

Definition at line 76 of file solver_2d.f90.

real*8, dimension(:,:,:), allocatable solver_2d::a_interface_xneg

Local speeds at the left of the x-interface.

Definition at line 47 of file solver_2d.f90.

real*8, dimension(:,:,:), allocatable solver_2d::a_interface_xpos

Local speeds at the right of the x-interface.

Definition at line 49 of file solver_2d.f90.

real*8, dimension(:,:,:), allocatable solver_2d::a_interface_yneg

Local speeds at the bottom of the y-interface.

Definition at line 51 of file solver_2d.f90.

real*8, dimension(:,:,:), allocatable solver_2d::a_interface_ypos

Local speeds at the top of the y-interface.

Definition at line 53 of file solver_2d.f90.

real*8, dimension(:), allocatable solver_2d::a_tilde

Explicit coeff. for the hyperbolic part for a single step of the R-K scheme.

Definition at line 84 of file solver_2d.f90.

real*8, dimension(:,:), allocatable solver_2d::a_tilde_ij

Butcher Tableau for the explicit part of the Runge-Kutta scheme.

Definition at line 74 of file solver_2d.f90.

real*8, dimension(:,:,:,:), allocatable solver_2d::divflux

Intermediate hyperbolic terms of the Runge-Kutta scheme.

Definition at line 93 of file solver_2d.f90.

real*8, dimension(:,:), allocatable solver_2d::divfluxj

Local Intermediate hyperbolic terms of the Runge-Kutta scheme.

Definition at line 103 of file solver_2d.f90.

real*8 solver_2d::dt

Time step.

Definition at line 67 of file solver_2d.f90.

real*8, dimension(:,:,:,:), allocatable solver_2d::expl_terms

Intermediate explicit terms of the Runge-Kutta scheme.

Definition at line 100 of file solver_2d.f90.

real*8, dimension(:,:), allocatable solver_2d::expl_terms_j

Local Intermediate explicit terms of the Runge-Kutta scheme.

Definition at line 107 of file solver_2d.f90.

real*8, dimension(:,:,:), allocatable solver_2d::h_interface_x

Semidiscrete numerical interface fluxes.

Definition at line 55 of file solver_2d.f90.

real*8, dimension(:,:,:), allocatable solver_2d::h_interface_y

Semidiscrete numerical interface fluxes.

Definition at line 57 of file solver_2d.f90.

integer solver_2d::i_rk

loop counter for the RK iteration

Definition at line 71 of file solver_2d.f90.

logical, dimension(:,:), allocatable solver_2d::mask11

Definition at line 69 of file solver_2d.f90.

logical, dimension(:,:), allocatable solver_2d::mask12

Definition at line 69 of file solver_2d.f90.

logical, dimension(:,:), allocatable solver_2d::mask21

Definition at line 69 of file solver_2d.f90.

logical, dimension(:,:), allocatable solver_2d::mask22

Definition at line 69 of file solver_2d.f90.

real*8, dimension(:,:,:,:), allocatable solver_2d::nh

Intermediate non-hyperbolic terms of the Runge-Kutta scheme.

Definition at line 95 of file solver_2d.f90.

real*8, dimension(:,:), allocatable solver_2d::nhj

Local Intermediate non-hyperbolic terms of the Runge-Kutta scheme.

Definition at line 105 of file solver_2d.f90.

logical solver_2d::normalize_f

Flag for the normalization of the array f in the implicit solution scheme.

Definition at line 115 of file solver_2d.f90.

logical solver_2d::normalize_q

Flag for the normalization of the array q in the implicit solution scheme.

Definition at line 112 of file solver_2d.f90.

real*8, dimension(:), allocatable solver_2d::omega

Coefficients for the implicit part of the Runge-Kutta scheme.

Definition at line 81 of file solver_2d.f90.

real*8, dimension(:), allocatable solver_2d::omega_tilde

Coefficients for the explicit part of the Runge-Kutta scheme.

Definition at line 79 of file solver_2d.f90.

logical solver_2d::opt_search_nl

Flag for the search of optimal step size in the implicit solution scheme.

Definition at line 118 of file solver_2d.f90.

real*8, dimension(:,:,:), allocatable solver_2d::q

Conservative variables.

Definition at line 33 of file solver_2d.f90.

real*8, dimension(:,:,:), allocatable solver_2d::q0

Conservative variables.

Definition at line 35 of file solver_2d.f90.

real*8, dimension(:,:,:), allocatable solver_2d::q_fv

Solution of the finite-volume semidiscrete cheme.

Definition at line 37 of file solver_2d.f90.

real*8, dimension(:,:,:), allocatable solver_2d::q_interfaceb

Reconstructed value at the bottom of the y-interface.

Definition at line 43 of file solver_2d.f90.

real*8, dimension(:,:,:), allocatable solver_2d::q_interfacel

Reconstructed value at the left of the x-interface.

Definition at line 39 of file solver_2d.f90.

real*8, dimension(:,:,:), allocatable solver_2d::q_interfacer

Reconstructed value at the right of the x-interface.

Definition at line 41 of file solver_2d.f90.

real*8, dimension(:,:,:), allocatable solver_2d::q_interfacet

Reconstructed value at the top of the y-interface.

Definition at line 45 of file solver_2d.f90.

real*8, dimension(:,:,:,:), allocatable solver_2d::q_rk

Intermediate solutions of the Runge-Kutta scheme.

Definition at line 91 of file solver_2d.f90.

real*8, dimension(:,:,:), allocatable solver_2d::qp

Physical variables ( $\alpha_1, p_1, p_2, \rho u, w, T$)

Definition at line 59 of file solver_2d.f90.

real*8, dimension(:,:,:), allocatable solver_2d::residual_term

Sum of all the terms of the equations except the transient term.

Definition at line 121 of file solver_2d.f90.

real*8, dimension(:,:,:,:), allocatable solver_2d::si_nh

Intermediate semi-implicit non-hyperbolic terms of the Runge-Kutta scheme.

Definition at line 97 of file solver_2d.f90.

real*8, dimension(:,:), allocatable solver_2d::si_nhj

Local Intermediate semi-impl non-hyperbolic terms of the Runge-Kutta scheme.

Definition at line 109 of file solver_2d.f90.

logical, dimension(:,:), allocatable solver_2d::solve_mask

Definition at line 64 of file solver_2d.f90.

logical, dimension(:,:), allocatable solver_2d::solve_mask0

Definition at line 64 of file solver_2d.f90.

real*8, dimension(:,:), allocatable solver_2d::source_xy

Array defining fraction of cells affected by source term.

Definition at line 62 of file solver_2d.f90.

real*8 solver_2d::t_imex1

Definition at line 124 of file solver_2d.f90.

real*8 solver_2d::t_imex2

Definition at line 124 of file solver_2d.f90.