SW_VAR_DENS_MODEL  0.9
Dept-averagedgas-particlemodel
Functions/Subroutines | Variables
constitutive_2d Module Reference

Constitutive equations. More...

Functions/Subroutines

subroutine init_problem_param
 Initialization of relaxation flags. More...
 
subroutine r_phys_var (r_qj, r_h, r_u, r_v, r_alphas, r_rho_m, r_T, r_alphal)
 Physical variables. More...
 
subroutine c_phys_var (c_qj, h, u, v, T, rho_m, red_grav, alphas)
 Physical variables. More...
 
subroutine mixt_var (qpj, r_Ri, r_rho_m, r_rho_c, r_red_grav)
 Physical variables. More...
 
subroutine qc_to_qp (qc, qp)
 Conservative to physical variables. More...
 
subroutine qp_to_qc (qp, B, qc)
 Physical to conservative variables. More...
 
subroutine qp_to_qp2 (qpj, Bj, qp2j)
 Additional Physical variables. More...
 
subroutine eval_local_speeds_x (qpj, vel_min, vel_max)
 Local Characteristic speeds x direction. More...
 
subroutine eval_local_speeds_y (qpj, vel_min, vel_max)
 Local Characteristic speeds y direction. More...
 
subroutine eval_fluxes (qcj, qpj, dir, flux)
 Hyperbolic Fluxes. More...
 
subroutine eval_nonhyperbolic_terms (c_qj, c_nh_term_impl, r_qj, r_nh_term_impl)
 Non-Hyperbolic terms. More...
 
subroutine eval_nh_semi_impl_terms (grav3_surf, qcj, nh_semi_impl_term)
 Non-Hyperbolic semi-implicit terms. More...
 
subroutine eval_expl_terms (Bprimej_x, Bprimej_y, source_xy, qpj, qcj, expl_term)
 Explicit Forces term. More...
 
subroutine eval_erosion_dep_term (qpj, dt, erosion_term, deposition_term)
 Erosion/Deposition term. More...
 
subroutine eval_topo_term (qpj, deposition_avg_term, erosion_avg_term, eqns_term, deposit_term)
 Topography modification related term. More...
 
subroutine eval_source_bdry (time, vect_x, vect_y, source_bdry)
 Internal boundary source fluxes. More...
 
real *8 function settling_velocity (diam, rhos, rhoc, i_solid)
 Settling velocity function. More...
 

Variables

logical, dimension(:), allocatable implicit_flag
 flag used for size of implicit non linear-system More...
 
logical entrainment_flag
 flag to activate air entrainment More...
 
real *8 grav
 gravitational acceleration More...
 
real *8 mu
 drag coefficients (Voellmy-Salm model) More...
 
real *8 xi
 
real *8 friction_factor
 drag coefficients (B&W model) More...
 
real *8 tau
 drag coefficients (plastic model) More...
 
real *8 t_env
 evironment temperature [K] More...
 
real *8 rad_coeff
 radiative coefficient More...
 
real *8 frict_coeff
 friction coefficient More...
 
real *8 t_ref
 reference temperature [K] More...
 
real *8 nu_ref
 reference kinematic viscosity [m2/s] More...
 
real *8 visc_par
 viscosity parameter [K-1] (b in Table 1 Costa & Macedonio, 2005) More...
 
real *8 emme
 velocity boundary layer fraction of total thickness More...
 
real *8 c_p
 specific heat [J kg-1 K-1] More...
 
real *8 atm_heat_transf_coeff
 atmospheric heat trasnfer coefficient [W m-2 K-1] (lambda in C&M, 2005) More...
 
real *8 exp_area_fract
 fractional area of the exposed inner core (f in C&M, 2005) More...
 
real *8, parameter sbconst = 5.67D-8
 Stephan-Boltzmann constant [W m-2 K-4]. More...
 
real *8 emissivity
 emissivity (eps in Costa & Macedonio, 2005) More...
 
real *8 enne
 thermal boundary layer fraction of total thickness More...
 
real *8 t_ground
 temperature of lava-ground interface More...
 
real *8 thermal_conductivity
 thermal conductivity [W m-1 K-1] (k in Costa & Macedonio, 2005) More...
 
real *8 alpha2
 1st param for yield strenght empirical relationship (O'Brian et al, 1993) More...
 
real *8 beta2
 2nd param for yield strenght empirical relationship (O'Brian et al, 1993) More...
 
real *8 alpha1_coeff
 ratio between reference value from input and computed values from eq. More...
 
real *8 beta1
 2nd param for fluid viscosity empirical relationship (O'Brian et al, 1993) More...
 
real *8 kappa
 Empirical resistance parameter (dimensionless, input parameter) More...
 
real *8 n_td
 Mannings roughness coefficient ( units: T L^(-1/3) ) More...
 
real *8 sp_heat_c
 Specific heat of carrier phase (gas or liquid) More...
 
real *8 rho_a_amb
 Ambient density of air ( units: kg m-3 ) More...
 
real *8 sp_heat_a
 Specific heat of air (units: J K-1 kg-1) More...
 
real *8 sp_gas_const_a
 Specific gas constant of air (units: J kg-1 K-1) More...
 
real *8 kin_visc_a
 Kinematic viscosity of air (units: m2 s-1) More...
 
real *8 kin_visc_l
 Kinematic viscosity of liquid (units: m2 s-1) More...
 
real *8 kin_visc_c
 Kinematic viscosity of carrier phase (units: m2 s-1) More...
 
real *8 t_ambient
 Temperature of ambient air (units: K) More...
 
real *8, dimension(:), allocatable rho_s
 Density of sediments ( units: kg m-3 ) More...
 
real *8, dimension(:), allocatable diam_s
 Diameter of sediments ( units: m ) More...
 
real *8, dimension(:), allocatable sp_heat_s
 Specific heat of solids (units: J K-1 kg-1) More...
 
logical settling_flag
 Flag to determine if sedimentation is active. More...
 
real *8 settling_vel
 Hindered settling velocity (units: m s-1 ) More...
 
real *8, dimension(:), allocatable erosion_coeff
 erosion model coefficient (units: m-1 ) More...
 
real *8 t_s_substrate
 temperature of solid substrate (units: K) More...
 
real *8 pres
 ambient pressure (units: Pa) More...
 
real *8 rho_l
 liquid density (units: kg m-3) More...
 
real *8 sp_heat_l
 Sepcific heat of liquid (units: J K-1 kg-1) More...
 

Detailed Description

Constitutive equations.

Function/Subroutine Documentation

subroutine constitutive_2d::c_phys_var ( complex*16, dimension(n_vars), intent(in)  c_qj,
complex*16, intent(out)  h,
complex*16, intent(out)  u,
complex*16, intent(out)  v,
complex*16, intent(out)  T,
complex*16, intent(out)  rho_m,
complex*16, intent(out)  red_grav,
complex*16, dimension(n_solid), intent(out)  alphas 
)

Physical variables.

This subroutine evaluates from the conservative local variables qj the local physical variables ( $h,u,v,T,\rho_m,red grav,\alpha_s $).

Parameters
[in]c_qjcomplex conservative variables
[out]hcomplex-value flow thickness
[out]ucomplex-value flow x-velocity
[out]vcomplex-value flow y-velocity
[out]Tcomplex-value flow temperature
[out]rho_mcomplex-value flow density
[out]red_gravcomplex-value flow density
[out]alphascomplex-value solid volume fractions
Author
Mattia de' Michieli Vitturi
Date
2019/12/13
Parameters
[out]hheight [m]
[out]uvelocity (x direction) [m s-1]
[out]vvelocity (y direction) [m s-1]
[out]ttemperature [K]
[out]rho_mmixture density [kg m-3]
[out]red_gravreduced gravity
[out]alphassediment volume fractions

Definition at line 388 of file constitutive_2d.f90.

Here is the caller graph for this function:

subroutine constitutive_2d::eval_erosion_dep_term ( real*8, dimension(n_vars+2), intent(in)  qpj,
real*8, intent(in)  dt,
real*8, dimension(n_solid), intent(out)  erosion_term,
real*8, dimension(n_solid), intent(out)  deposition_term 
)

Erosion/Deposition term.

This subroutine evaluates the deposition term.

Date
03/010/2018
Parameters
[in]qpjlocal physical variables
[in]dttime step
[out]erosion_termerosion term for each solid phase
[out]dep_termdeposition term for each solid phase
Author
Mattia de' Michieli Vitturi
Parameters
[in]qpjphysical variables
[out]erosion_termerosion term
[out]deposition_termdeposition term

Definition at line 1668 of file constitutive_2d.f90.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine constitutive_2d::eval_expl_terms ( real*8, intent(in)  Bprimej_x,
real*8, intent(in)  Bprimej_y,
real*8, intent(in)  source_xy,
real*8, dimension(n_vars+2), intent(in)  qpj,
real*8, dimension(n_vars), intent(in)  qcj,
real*8, dimension(n_eqns), intent(out)  expl_term 
)

Explicit Forces term.

This subroutine evaluates the non-hyperbolic terms to be treated explicitely in the DIRK numerical scheme (e.g. gravity,source of mass). The sign of the terms is taken with the terms on the left-hand side of the equations.

Date
2019/12/13
Parameters
[in]B_primej_xlocal x-slope
[in]B_primej_ylocal y_slope
[in]source_xylocal source
[in]qpjphysical variables
[in]qcjconservative variables
[out]expl_termexplicit term
Author
Mattia de' Michieli Vitturi
Parameters
[in]qpjlocal physical variables
[in]qcjlocal conservative variables
[out]expl_termlocal explicit forces

Definition at line 1604 of file constitutive_2d.f90.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine constitutive_2d::eval_fluxes ( real*8, dimension(n_vars), intent(in)  qcj,
real*8, dimension(n_vars+2), intent(in)  qpj,
integer, intent(in)  dir,
real*8, dimension(n_eqns), intent(out)  flux 
)

Hyperbolic Fluxes.

This subroutine evaluates the numerical fluxes given the conservative variables qcj and physical variables qpj.

Date
01/06/2012
Parameters
[in]qcjreal local conservative variables
[in]qpjreal local physical variables
[in]dirdirection of the flux (1=x,2=y)
[out]fluxreal fluxes
Author
Mattia de' Michieli Vitturi

Definition at line 1077 of file constitutive_2d.f90.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine constitutive_2d::eval_local_speeds_x ( real*8, dimension(n_vars+2), intent(in)  qpj,
real*8, dimension(n_vars), intent(out)  vel_min,
real*8, dimension(n_vars), intent(out)  vel_max 
)

Local Characteristic speeds x direction.

This subroutine computes from the physical variable evaluates the largest positive and negative characteristic speed in the x-direction.

Parameters
[in]qpjarray of local physical variables
[out]vel_minminimum x-velocity
[out]vel_maxmaximum x-velocity
Author
Mattia de' Michieli Vitturi
Date
05/12/2017

Definition at line 973 of file constitutive_2d.f90.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine constitutive_2d::eval_local_speeds_y ( real*8, dimension(n_vars+2), intent(in)  qpj,
real*8, dimension(n_vars), intent(out)  vel_min,
real*8, dimension(n_vars), intent(out)  vel_max 
)

Local Characteristic speeds y direction.

This subroutine computes from the physical variable evaluates the largest positive and negative characteristic speed in the y-direction.

Parameters
[in]qpjarray of local physical variables
[out]vel_minminimum y-velocity
[out]vel_maxmaximum y-velocity
Author
Mattia de' Michieli Vitturi
Date
05/12/2017

Definition at line 1024 of file constitutive_2d.f90.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine constitutive_2d::eval_nh_semi_impl_terms ( real*8, intent(in)  grav3_surf,
real*8, dimension(n_vars), intent(in)  qcj,
real*8, dimension(n_eqns), intent(out)  nh_semi_impl_term 
)

Non-Hyperbolic semi-implicit terms.

This subroutine evaluates the non-hyperbolic terms that are solved semi-implicitely by the solver. For example, any discontinuous term that appears in the friction terms.

Date
20/01/2018
Parameters
[in]grav3_surfgravity correction
[in]qcjreal conservative variables
[out]nh_semi_impl_termreal non-hyperbolic terms
Author
Mattia de' Michieli Vitturi

Yield strenght (units: kg m−1 s−2)

Yield slope component of total friction (dimensionless)

Definition at line 1478 of file constitutive_2d.f90.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine constitutive_2d::eval_nonhyperbolic_terms ( complex*16, dimension(n_vars), intent(in), optional  c_qj,
complex*16, dimension(n_eqns), intent(out), optional  c_nh_term_impl,
real*8, dimension(n_vars), intent(in), optional  r_qj,
real*8, dimension(n_eqns), intent(out), optional  r_nh_term_impl 
)

Non-Hyperbolic terms.

This subroutine evaluates the non-hyperbolic terms (relaxation terms and forces) of the system of equations, both for real or complex inputs. These terms are treated implicitely in the DIRK numerical scheme.

Date
01/06/2012
Parameters
[in]c_qjcomplex conservative variables
[in]r_qjreal conservative variables
[out]c_nh_term_implcomplex non-hyperbolic terms
[out]r_nh_term_implreal non-hyperbolic terms
Author
Mattia de' Michieli Vitturi

Temperature in C

1st param for fluid viscosity empirical relationship (O'Brian et al, 1993)

Fluid dynamic viscosity (units: kg m-1 s-1 )

Total friction (dimensionless)

Viscous slope component of total Friction (dimensionless)

Turbulent dispersive slope component of total friction (dimensionless)

Definition at line 1208 of file constitutive_2d.f90.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine constitutive_2d::eval_source_bdry ( real*8, intent(in)  time,
real*8, intent(in)  vect_x,
real*8, intent(in)  vect_y,
real*8, dimension(n_vars), intent(out)  source_bdry 
)

Internal boundary source fluxes.

This subroutine evaluates the source terms at the interfaces when an internal radial source is present, as for a base surge. The terms are applied as boundary conditions, and thus they have the units of the physical variable qp

Date
2019/12/01
Parameters
[in]timetime
[in]vect_xunit vector velocity x-component
[in]vect_yunit vector velocity y-component
[out]source_bdrysource terms
Author
Mattia de' Michieli Vitturi

Definition at line 1888 of file constitutive_2d.f90.

Here is the caller graph for this function:

subroutine constitutive_2d::eval_topo_term ( real*8, dimension(n_vars+2), intent(in)  qpj,
real*8, dimension(n_solid), intent(in)  deposition_avg_term,
real*8, dimension(n_solid), intent(in)  erosion_avg_term,
real*8, dimension(n_eqns), intent(out)  eqns_term,
real*8, dimension(n_solid), intent(out)  deposit_term 
)

Topography modification related term.

This subroutine evaluates the deposition term.

Date
2019/11/08
Parameters
[in]qpjphysical variables
[in]deposition_avg_termaveraged deposition terms
[in]erosion_avg_termaveraged deposition terms
[out]eqns_termsource terms for cons equations
[out]deposit_termdeposition rates for solids
Author
Mattia de' Michieli Vitturi
Parameters
[in]qpjphysical variables
[in]deposition_avg_termdeposition term
[in]erosion_avg_termerosion term

Definition at line 1774 of file constitutive_2d.f90.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine constitutive_2d::init_problem_param ( )

Initialization of relaxation flags.

This subroutine set the number and the flags of the non-hyperbolic terms.

Date
07/09/2012

Definition at line 164 of file constitutive_2d.f90.

Here is the caller graph for this function:

subroutine constitutive_2d::mixt_var ( real*8, dimension(n_vars+2), intent(in)  qpj,
real*8, intent(out)  r_Ri,
real*8, intent(out)  r_rho_m,
real*8, intent(out)  r_rho_c,
real*8, intent(out)  r_red_grav 
)

Physical variables.

This subroutine evaluates from the physical real-value local variables qpj, all the (real-valued ) variables that define the physical state and that are needed to compute the explicit equations terms.

Parameters
[in]qpjreal-valued physical variables
[out]r_Rireal-valued Richardson number
[out]r_rho_mreal-valued mixture density
[out]r_rho_creal-valued carrier phase density
[out]r_red_gravreal-valued reduced gravity
Author
Mattia de' Michieli Vitturi
Date
10/10/2019
Parameters
[in]qpjreal-value physical variables
[out]r_rireal-value Richardson number
[out]r_rho_mreal-value mixture density [kg/m3]
[out]r_rho_creal-value carrier phase density [kg/m3]
[out]r_red_gravreal-value reduced gravity

Definition at line 569 of file constitutive_2d.f90.

Here is the caller graph for this function:

subroutine constitutive_2d::qc_to_qp ( real*8, dimension(n_vars), intent(in)  qc,
real*8, dimension(n_vars+2), intent(out)  qp 
)

Conservative to physical variables.

This subroutine evaluates from the conservative variables qc the array of physical variables qp:

  • qp(1) = $ h $
  • qp(2) = $ hu $
  • qp(3) = $ hv $
  • qp(4) = $ T $
  • qp(5:4+n_solid) = $ alphas(1:n_solid) $
  • qp(n_vars) = $ alphal $
  • qp(n_vars+1) = $ u $
  • qp(n_vars+2) = $ v $

The physical variables are those used for the linear reconstruction at the cell interfaces. Limiters are applied to the reconstructed slopes.

Parameters
[in]qclocal conservative variables
[out]qplocal physical variables
Date
2019/11/11
Author
Mattia de' Michieli Vitturi

Definition at line 678 of file constitutive_2d.f90.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine constitutive_2d::qp_to_qc ( real*8, dimension(n_vars+2), intent(in)  qp,
real*8, intent(in)  B,
real*8, dimension(n_vars), intent(out)  qc 
)

Physical to conservative variables.

This subroutine evaluates the conservative real_value variables qc from the array of real_valued physical variables qp:

  • qp(1) = $ h $
  • qp(2) = $ h*u $
  • qp(3) = $ h*v $
  • qp(4) = $ T $
  • qp(5:4+n_s) = $ alphas(1:n_s) $
  • qp(n_vars) = $ alphal $
  • qp(n_vars+1) = $ u $
  • qp(n_vars+2) = $ v $
Parameters
[in]qpphysical variables
[in]Blocal topography
[out]qcconservative variables
Date
2019/11/18
Author
Mattia de' Michieli Vitturi

Definition at line 752 of file constitutive_2d.f90.

Here is the caller graph for this function:

subroutine constitutive_2d::qp_to_qp2 ( real*8, dimension(n_vars), intent(in)  qpj,
real*8, intent(in)  Bj,
real*8, dimension(3), intent(out)  qp2j 
)

Additional Physical variables.

This subroutine evaluates from the physical local variables qpj, the two additional local variables qp2j = (h+B,u,v).

Parameters
[in]qpjreal-valued physical variables
[in]Bjreal-valued local topography
[out]qp2jreal-valued physical variables
Author
Mattia de' Michieli Vitturi
Date
10/10/2019

Definition at line 934 of file constitutive_2d.f90.

Here is the caller graph for this function:

subroutine constitutive_2d::r_phys_var ( real*8, dimension(n_vars), intent(in)  r_qj,
real*8, intent(out)  r_h,
real*8, intent(out)  r_u,
real*8, intent(out)  r_v,
real*8, dimension(n_solid), intent(out)  r_alphas,
real*8, intent(out)  r_rho_m,
real*8, intent(out)  r_T,
real*8, intent(out)  r_alphal 
)

Physical variables.

This subroutine evaluates from the conservative local variables qj the local physical variables ( $h,u,v,\alpha_s,\rho_m,T,\alpha_l $).

Parameters
[in]r_qjreal conservative variables
[out]r_hreal-value flow thickness
[out]r_ureal-value flow x-velocity
[out]r_vreal-value flow y-velocity
[out]r_alphasreal-value solid volume fractions
[out]r_rho_mreal-value flow density
[out]r_Treal-value flow temperature
[out]r_alphalreal-value liquid volume fraction
Author
Mattia de' Michieli Vitturi
Date
2019/12/13
Parameters
[in]r_qjreal-value conservative var
[out]r_hreal-value flow thickness
[out]r_ureal-value x-velocity
[out]r_vreal-value y-velocity
[out]r_alphasreal-value solid volume fractions
[out]r_rho_mreal-value mixture density
[out]r_treal-value temperature
[out]r_alphalreal-value liquid volume fraction

Definition at line 206 of file constitutive_2d.f90.

Here is the caller graph for this function:

real*8 function constitutive_2d::settling_velocity ( real*8, intent(in)  diam,
real*8, intent(in)  rhos,
real*8, intent(in)  rhoc,
integer, intent(in)  i_solid 
)

Settling velocity function.

This subroutine compute the settling velocity of the particles, as a function of diameter, density.

Date
2019/11/11
Parameters
[in]diamparticle diameter
[in]rhosparticle density
[in]rhoaatmospheric density [in] i_solid particle class index
Author
Mattia de' Michieli Vitturi
Parameters
[in]diamparticle diameter [m]
[in]rhosparticle density [kg/m3]
[in]rhoccarrier phase density [kg/m3]
[in]i_solidparticle class index

Definition at line 1991 of file constitutive_2d.f90.

Here is the caller graph for this function:

Variable Documentation

real*8 constitutive_2d::alpha1_coeff

ratio between reference value from input and computed values from eq.

Definition at line 86 of file constitutive_2d.f90.

real*8 constitutive_2d::alpha2

1st param for yield strenght empirical relationship (O'Brian et al, 1993)

Definition at line 80 of file constitutive_2d.f90.

real*8 constitutive_2d::atm_heat_transf_coeff

atmospheric heat trasnfer coefficient [W m-2 K-1] (lambda in C&M, 2005)

Definition at line 57 of file constitutive_2d.f90.

real*8 constitutive_2d::beta1

2nd param for fluid viscosity empirical relationship (O'Brian et al, 1993)

Definition at line 89 of file constitutive_2d.f90.

real*8 constitutive_2d::beta2

2nd param for yield strenght empirical relationship (O'Brian et al, 1993)

Definition at line 83 of file constitutive_2d.f90.

real*8 constitutive_2d::c_p

specific heat [J kg-1 K-1]

Definition at line 54 of file constitutive_2d.f90.

real*8, dimension(:), allocatable constitutive_2d::diam_s

Diameter of sediments ( units: m )

Definition at line 127 of file constitutive_2d.f90.

real*8 constitutive_2d::emissivity

emissivity (eps in Costa & Macedonio, 2005)

Definition at line 66 of file constitutive_2d.f90.

real*8 constitutive_2d::emme

velocity boundary layer fraction of total thickness

Definition at line 51 of file constitutive_2d.f90.

real*8 constitutive_2d::enne

thermal boundary layer fraction of total thickness

Definition at line 69 of file constitutive_2d.f90.

logical constitutive_2d::entrainment_flag

flag to activate air entrainment

Definition at line 17 of file constitutive_2d.f90.

real*8, dimension(:), allocatable constitutive_2d::erosion_coeff

erosion model coefficient (units: m-1 )

Definition at line 139 of file constitutive_2d.f90.

real*8 constitutive_2d::exp_area_fract

fractional area of the exposed inner core (f in C&M, 2005)

Definition at line 60 of file constitutive_2d.f90.

real*8 constitutive_2d::frict_coeff

friction coefficient

Definition at line 39 of file constitutive_2d.f90.

real*8 constitutive_2d::friction_factor

drag coefficients (B&W model)

Definition at line 27 of file constitutive_2d.f90.

real*8 constitutive_2d::grav

gravitational acceleration

Definition at line 20 of file constitutive_2d.f90.

logical, dimension(:), allocatable constitutive_2d::implicit_flag

flag used for size of implicit non linear-system

Definition at line 14 of file constitutive_2d.f90.

real*8 constitutive_2d::kappa

Empirical resistance parameter (dimensionless, input parameter)

Definition at line 92 of file constitutive_2d.f90.

real*8 constitutive_2d::kin_visc_a

Kinematic viscosity of air (units: m2 s-1)

Definition at line 112 of file constitutive_2d.f90.

real*8 constitutive_2d::kin_visc_c

Kinematic viscosity of carrier phase (units: m2 s-1)

Definition at line 118 of file constitutive_2d.f90.

real*8 constitutive_2d::kin_visc_l

Kinematic viscosity of liquid (units: m2 s-1)

Definition at line 115 of file constitutive_2d.f90.

real*8 constitutive_2d::mu

drag coefficients (Voellmy-Salm model)

Definition at line 23 of file constitutive_2d.f90.

real*8 constitutive_2d::n_td

Mannings roughness coefficient ( units: T L^(-1/3) )

Definition at line 95 of file constitutive_2d.f90.

real*8 constitutive_2d::nu_ref

reference kinematic viscosity [m2/s]

Definition at line 45 of file constitutive_2d.f90.

real*8 constitutive_2d::pres

ambient pressure (units: Pa)

Definition at line 145 of file constitutive_2d.f90.

real*8 constitutive_2d::rad_coeff

radiative coefficient

Definition at line 36 of file constitutive_2d.f90.

real*8 constitutive_2d::rho_a_amb

Ambient density of air ( units: kg m-3 )

Definition at line 103 of file constitutive_2d.f90.

real*8 constitutive_2d::rho_l

liquid density (units: kg m-3)

Definition at line 148 of file constitutive_2d.f90.

real*8, dimension(:), allocatable constitutive_2d::rho_s

Density of sediments ( units: kg m-3 )

Definition at line 124 of file constitutive_2d.f90.

real*8, parameter constitutive_2d::sbconst = 5.67D-8

Stephan-Boltzmann constant [W m-2 K-4].

Definition at line 63 of file constitutive_2d.f90.

logical constitutive_2d::settling_flag

Flag to determine if sedimentation is active.

Definition at line 133 of file constitutive_2d.f90.

real*8 constitutive_2d::settling_vel

Hindered settling velocity (units: m s-1 )

Definition at line 136 of file constitutive_2d.f90.

real*8 constitutive_2d::sp_gas_const_a

Specific gas constant of air (units: J kg-1 K-1)

Definition at line 109 of file constitutive_2d.f90.

real*8 constitutive_2d::sp_heat_a

Specific heat of air (units: J K-1 kg-1)

Definition at line 106 of file constitutive_2d.f90.

real*8 constitutive_2d::sp_heat_c

Specific heat of carrier phase (gas or liquid)

Definition at line 100 of file constitutive_2d.f90.

real*8 constitutive_2d::sp_heat_l

Sepcific heat of liquid (units: J K-1 kg-1)

Definition at line 151 of file constitutive_2d.f90.

real*8, dimension(:), allocatable constitutive_2d::sp_heat_s

Specific heat of solids (units: J K-1 kg-1)

Definition at line 130 of file constitutive_2d.f90.

real*8 constitutive_2d::t_ambient

Temperature of ambient air (units: K)

Definition at line 121 of file constitutive_2d.f90.

real*8 constitutive_2d::t_env

evironment temperature [K]

Definition at line 33 of file constitutive_2d.f90.

real*8 constitutive_2d::t_ground

temperature of lava-ground interface

Definition at line 72 of file constitutive_2d.f90.

real*8 constitutive_2d::t_ref

reference temperature [K]

Definition at line 42 of file constitutive_2d.f90.

real*8 constitutive_2d::t_s_substrate

temperature of solid substrate (units: K)

Definition at line 142 of file constitutive_2d.f90.

real*8 constitutive_2d::tau

drag coefficients (plastic model)

Definition at line 30 of file constitutive_2d.f90.

real*8 constitutive_2d::thermal_conductivity

thermal conductivity [W m-1 K-1] (k in Costa & Macedonio, 2005)

Definition at line 75 of file constitutive_2d.f90.

real*8 constitutive_2d::visc_par

viscosity parameter [K-1] (b in Table 1 Costa & Macedonio, 2005)

Definition at line 48 of file constitutive_2d.f90.

real*8 constitutive_2d::xi

Definition at line 24 of file constitutive_2d.f90.