LAHARS-MODEL  0.1
templategithubproject
Functions/Subroutines | Variables
constitutive_2d Module Reference

Constitutive equations. More...

Functions/Subroutines

subroutine init_problem_param
 Initialization of relaxation flags. More...
 
subroutine phys_var (r_qj, c_qj)
 Physical variables. More...
 
subroutine eval_local_speeds_x (qj, vel_min, vel_max)
 Local Characteristic speeds x direction. More...
 
subroutine eval_local_speeds_y (qj, vel_min, vel_max)
 Local Characteristic speeds y direction. More...
 
subroutine qc_to_qp (qc, B, qp)
 Conservative to physical variables. More...
 
subroutine qp_to_qc (qp, B, qc)
 Physical to conservative variables. More...
 
subroutine qc_to_qc2 (qc, B, qc2)
 Conservative to alternative conservative variables. More...
 
subroutine qc2_to_qc (qc2, B, qc)
 Reconstructed to conservative variables. More...
 
subroutine eval_fluxes (c_qj, r_qj, c_flux, r_flux, dir)
 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, c_qj, c_nh_semi_impl_term, r_qj, r_nh_semi_impl_term)
 Non-Hyperbolic semi-implicit terms. More...
 
subroutine eval_expl_terms (Bprimej_x, Bprimej_y, source_xy, qj, expl_term)
 Explicit Forces term. More...
 
subroutine eval_erosion_dep_term (qj, erosion_term, deposition_term)
 Deposition term. More...
 
subroutine eval_topo_term (qj, deposition_avg_term, erosion_avg_term, eqns_term, topo_term)
 Topography modification related term. More...
 

Variables

logical, dimension(:), allocatable implicit_flag
 
character(len=20) phase1_name
 
character(len=20) phase2_name
 
complex *16 h
 height [m] More...
 
complex *16 u
 velocity (x direction) More...
 
complex *16 v
 velocity (y direction) More...
 
complex *16 t
 temperature More...
 
complex *16 alphas
 sediment volume fraction More...
 
complex *16 rho_m
 mixture density More...
 
real *8 grav
 gravitational acceleration More...
 
real *8 mu
 drag coefficients (Voellmy-Salm model) More...
 
real *8 xi
 
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 rho
 fluid density [kg/m3] 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 parameter for yield strenght empirical relationship (O'Brian et al, 1993) More...
 
real *8 beta2
 2nd parameter for yield strenght empirical relationship (O'Brian et al, 1993) More...
 
real *8 alpha1
 1st parameter for fluid viscosity empirical relationship (O'Brian et al, 1993) More...
 
real *8 beta1
 2nd parameter for fluid viscosity empirical relationship (O'Brian et al, 1993) More...
 
real *8 kappa
 Empirical resistance parameter. More...
 
real *8 n_td
 Mannings roughness coefficient ( units: T L^(-1/3) ) More...
 
real *8 rho_w
 Specific weight of water. More...
 
real *8 rho_s
 Specific weight of sediments. More...
 
real *8 settling_vel
 hindered settling More...
 
real *8 erosion_coeff
 erosion model coefficient More...
 

Detailed Description

Constitutive equations.

Function/Subroutine Documentation

subroutine constitutive_2d::eval_erosion_dep_term ( real*8, dimension(n_eqns), intent(in)  qj,
real*8, intent(out)  erosion_term,
real*8, intent(out)  deposition_term 
)

Deposition term.

This subroutine evaluates the deposition term.

Date
03/010/2018
Parameters
[in]Bjtopography
[in]qjconservative variables
[out]dep_termexplicit term
[in]qjconservative variables
[out]erosion_termerosion term
[out]deposition_termdeposition term

settling speed of the single grain

Richardson-Zaki exponent [Richardson and Zaki, 1954]

settling particle Reynold number

kinematic viscosity of the fluid

grain diameter

dimensionless grain parameter

Definition at line 1283 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_eqns), intent(in)  qj,
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
01/06/2012
Parameters
[in]qjconservative variables
[out]expl_termexplicit term
[in]qjconservative variables
[out]expl_termexplicit forces

Definition at line 1200 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 ( complex*16, dimension(n_vars), intent(in), optional  c_qj,
real*8, dimension(n_vars), intent(in), optional  r_qj,
complex*16, dimension(n_eqns), intent(out), optional  c_flux,
real*8, dimension(n_eqns), intent(out), optional  r_flux,
integer, intent(in)  dir 
)

Hyperbolic Fluxes.

This subroutine evaluates the numerical fluxes given the conservative variables qj, accordingly to the equations for the single temperature model introduced in Romenki et al. 2010.

Date
01/06/2012
Parameters
[in]c_qjcomplex conservative variables
[in]r_qjreal conservative variables
[out]c_fluxcomplex analytical fluxes
[out]r_fluxreal analytical fluxes

Definition at line 563 of file constitutive_2d.f90.

Here is the caller graph for this function:

subroutine constitutive_2d::eval_local_speeds_x ( real*8, dimension(n_vars), intent(in)  qj,
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 desingularize the velocities and then evaluates the largest positive and negative characteristic speed in the x-direction.

Parameters
[in]qjarray of conservative variables
[in]Bjtopography at the cell center
[out]vel_minminimum x-velocity
[out]vel_maxmaximum x-velocity
Author
Mattia de' Michieli Vitturi
Date
05/12/2017

Definition at line 262 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), intent(in)  qj,
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 desingularize the velocities and then evaluates the largest positive and negative characteristic speed in the y-direction.

Parameters
[in]qjarray of conservative variables
[in]Bjtopography at the cell center
[out]vel_minminimum y-velocity
[out]vel_maxmaximum y-velocity
Author
Mattia de' Michieli Vitturi
Date
05/12/2017

Definition at line 291 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,
complex*16, dimension(n_vars), intent(in), optional  c_qj,
complex*16, dimension(n_eqns), intent(out), optional  c_nh_semi_impl_term,
real*8, dimension(n_vars), intent(in), optional  r_qj,
real*8, dimension(n_eqns), intent(out), optional  r_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]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

Yield strenght

Yield slope component of total friction

Definition at line 1048 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

Fluid viscosity

Total friction

Viscous slope component of total Friction

Turbulent dispersive slope component of total friction

Definition at line 725 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_topo_term ( real*8, dimension(n_eqns), intent(in)  qj,
real*8, intent(in)  deposition_avg_term,
real*8, intent(in)  erosion_avg_term,
real*8, dimension(n_eqns), intent(out)  eqns_term,
real*8, intent(out)  topo_term 
)

Topography modification related term.

This subroutine evaluates the deposition term.

Date
03/010/2018
Parameters
[in]Bjtopography
[in]qjconservative variables
[in]deposition_avg_termaveraged deposition terms
[in]erosion_avg_termaveraged deposition terms
[out]topo_termexplicit term
[in]qjconservative variables
[in]deposition_avg_termdeposition term
[in]erosion_avg_termerosion term

Definition at line 1371 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 128 of file constitutive_2d.f90.

Here is the caller graph for this function:

subroutine constitutive_2d::phys_var ( real*8, dimension(n_vars), intent(in), optional  r_qj,
complex*16, dimension(n_vars), intent(in), optional  c_qj 
)

Physical variables.

This subroutine evaluates from the conservative local variables qj the local physical variables ( $h+B, u, v, xs , T $).

Parameters
[in]r_qjreal conservative variables
[in]c_qjcomplex conservative variables
Author
Mattia de' Michieli Vitturi
Date
15/08/2011

Definition at line 164 of file constitutive_2d.f90.

Here is the caller graph for this function:

subroutine constitutive_2d::qc2_to_qc ( real*8, dimension(n_vars), intent(in)  qc2,
real*8, intent(in)  B,
real*8, dimension(n_vars), intent(out)  qc 
)

Reconstructed to conservative variables.

This subroutine evaluates the conservative variables qc from the array of alternative conser variables qc2, reconstructed at the interfaces:

  • qc2(1) = $ h + B $
  • qc2(2) = $ rho_m hu $
  • qc2(3) = $ rho_m hv $
  • qc2(4) = $ h \cdot alpha_s $
  • qc2(5) = $ h \cdot T $
Parameters
[in]qc2alternative conservative variables
[out]qcconservative variables
Date
02/04/2019

Definition at line 480 of file constitutive_2d.f90.

Here is the call graph for this function:

Here is the caller graph for this function:

subroutine constitutive_2d::qc_to_qc2 ( real*8, dimension(n_vars), intent(in)  qc,
real*8, intent(in)  B,
real*8, dimension(n_vars), intent(out)  qc2 
)

Conservative to alternative conservative variables.

This subroutine evaluates from the conservative variables qc the array of alternative conservative variables qc2:

  • qc2(1) = $ h+B $
  • qc2(2) = $ rho_m hu $
  • qc2(3) = $ rho_m hv $
  • qc2(4) = $ h alpha_s $
  • qc2(5) = $ hT $

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

Parameters
[in]qcconservative variables
[out]qc2alternative conservative variables
Date
02/04/2019

Definition at line 451 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, intent(in)  B,
real*8, dimension(n_vars), intent(out)  qp 
)

Conservative to physical variables.

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

  • qp(1) = $ h+B $
  • qp(2) = $ u $
  • qp(3) = $ v $
  • qp(4) = $ xs $
  • qp(5) = $ T $

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

Parameters
[in]qcconservative variables
[out]qpphysical variables
Date
15/08/2011

Definition at line 323 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), 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 variables qc from the array of physical variables qp:

  • qp(1) = $ h + B $
  • qp(2) = $ rho_m*u $
  • qp(3) = $ rho_m*v $
  • qp(4) = $ alphas $
  • qp(5) = $ T $
Parameters
[in]qpphysical variables
[out]qcconservative variables
Date
15/08/2011

Definition at line 372 of file constitutive_2d.f90.

Here is the caller graph for this function:

Variable Documentation

real*8 constitutive_2d::alpha1

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

Definition at line 93 of file constitutive_2d.f90.

real*8 constitutive_2d::alpha2

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

Definition at line 87 of file constitutive_2d.f90.

complex*16 constitutive_2d::alphas

sediment volume fraction

Definition at line 23 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 64 of file constitutive_2d.f90.

real*8 constitutive_2d::beta1

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

Definition at line 96 of file constitutive_2d.f90.

real*8 constitutive_2d::beta2

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

Definition at line 90 of file constitutive_2d.f90.

real*8 constitutive_2d::c_p

specific heat [J kg-1 K-1]

Definition at line 61 of file constitutive_2d.f90.

real*8 constitutive_2d::emissivity

emissivity (eps in Costa & Macedonio, 2005)

Definition at line 73 of file constitutive_2d.f90.

real*8 constitutive_2d::emme

velocity boundary layer fraction of total thickness

Definition at line 58 of file constitutive_2d.f90.

real*8 constitutive_2d::enne

thermal boundary layer fraction of total thickness

Definition at line 76 of file constitutive_2d.f90.

real*8 constitutive_2d::erosion_coeff

erosion model coefficient

Definition at line 114 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 67 of file constitutive_2d.f90.

real*8 constitutive_2d::frict_coeff

friction coefficient

Definition at line 43 of file constitutive_2d.f90.

real*8 constitutive_2d::grav

gravitational acceleration

Definition at line 27 of file constitutive_2d.f90.

complex*16 constitutive_2d::h

height [m]

Definition at line 19 of file constitutive_2d.f90.

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

Definition at line 14 of file constitutive_2d.f90.

real*8 constitutive_2d::kappa

Empirical resistance parameter.

Definition at line 99 of file constitutive_2d.f90.

real*8 constitutive_2d::mu

drag coefficients (Voellmy-Salm model)

Definition at line 30 of file constitutive_2d.f90.

real*8 constitutive_2d::n_td

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

Definition at line 102 of file constitutive_2d.f90.

real*8 constitutive_2d::nu_ref

reference kinematic viscosity [m2/s]

Definition at line 52 of file constitutive_2d.f90.

character(len=20) constitutive_2d::phase1_name

Definition at line 16 of file constitutive_2d.f90.

character(len=20) constitutive_2d::phase2_name

Definition at line 17 of file constitutive_2d.f90.

real*8 constitutive_2d::rad_coeff

radiative coefficient

Definition at line 40 of file constitutive_2d.f90.

real*8 constitutive_2d::rho

fluid density [kg/m3]

Definition at line 46 of file constitutive_2d.f90.

complex*16 constitutive_2d::rho_m

mixture density

Definition at line 24 of file constitutive_2d.f90.

real*8 constitutive_2d::rho_s

Specific weight of sediments.

Definition at line 108 of file constitutive_2d.f90.

real*8 constitutive_2d::rho_w

Specific weight of water.

Definition at line 105 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 70 of file constitutive_2d.f90.

real*8 constitutive_2d::settling_vel

hindered settling

Definition at line 111 of file constitutive_2d.f90.

complex*16 constitutive_2d::t

temperature

Definition at line 22 of file constitutive_2d.f90.

real*8 constitutive_2d::t_env

evironment temperature [K]

Definition at line 37 of file constitutive_2d.f90.

real*8 constitutive_2d::t_ground

temperature of lava-ground interface

Definition at line 79 of file constitutive_2d.f90.

real*8 constitutive_2d::t_ref

reference temperature [K]

Definition at line 49 of file constitutive_2d.f90.

real*8 constitutive_2d::tau

drag coefficients (plastic model)

Definition at line 34 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 82 of file constitutive_2d.f90.

complex*16 constitutive_2d::u

velocity (x direction)

Definition at line 20 of file constitutive_2d.f90.

complex*16 constitutive_2d::v

velocity (y direction)

Definition at line 21 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 55 of file constitutive_2d.f90.

real*8 constitutive_2d::xi

Definition at line 31 of file constitutive_2d.f90.