exo_k.atm_evolution.convection
@author: jeremy leconte
Module Contents
- exo_k.atm_evolution.convection.dry_convective_adjustment_numba(timestep, Nlay, t_lay, exner, dmass, tracer_array, Mgas, verbose=False)[source]
Computes the heating rates needed to adjust unstable regions of a given atmosphere to a convectively neutral T profile on a given timestep.
- Parameters:
timestep (
float
) – Duration of the adjustment. If given in seconds, H=dT/timestep is in K/s. In the current model, timestep is in second over cp. This ensures that H=dT/timestep is in W/kg.Nlay (
int
) – Number of atmospheric layerst_lay (
array
,np.ndarray
) – Temperatures of the atmospheric layers (K)exner (
array
,np.ndarray
) –Exner function computed at the layer centers ((p/psurf)**rcp)
\[\Pi=(p / p_{s})^{R/c_p}\]dmass (
array
,np.ndarray
) – Mass of gas in each layer (kg/m^2)
- Returns:
- array
Heating rate in each atmospheric layer (W/kg).
- exo_k.atm_evolution.convection.convective_acceleration_numba(timestep, Nlay, H_rad, rad_layers, tau_rad, tau_rads, dmass, convective_acceleration_mode=0, convective_acceleration_factor=1.0, verbose=False)[source]
Computes a heating rate for the whole convective region to accelerate convergence
- Parameters:
timestep (
float
) – Duration of the adjustment. If given in seconds, H=dT/timestep is in K/s. In the current model, timestep is in second over cp. This ensures that H=dT/timestep is in W/kg.Nlay (
int
) – Number of atmospheric layersH_rad (
array
,np.ndarray
) – Radiative heating raterad_layers (
array
,np.ndarray
ofbool
) – Elements in the array are true if layer is purely radiativetau_rad (
array
,np.ndarray
) – Baseline radiative timescale for the atmosphere. (e.g. the min of tau_rads) Should use the same units as timestep.tau_rads (
array
,np.ndarray
) – Radiative timescale for each layer. Should use the same units as timestep.dmass (
array
,np.ndarray
) – Mass of gas in each layer (kg/m^2)convective_acceleration_mode (
int
) – convective_acceleration_mode=0 is the default mode =1 emulates an old (and bad) behavior.convective_acceleration_factor (
int
) – Multiplying factor for convection acceleration
- Returns:
- array
Heating rate in each atmospheric layer (W/kg).
- exo_k.atm_evolution.convection.turbulent_diffusion_numba(timestep, Nlay, p_lay, p_lev, dmass, t_lay, exner, t_lay_ov_mu, g, Kzz, tracer_array, cp, mix_potential_temp=False, verbose=False)[source]
Solves turbulent diffusion equation:
\[\rho frac{\partial q}{\partial t} = \frac{\partial F_{diff}}{\partial z}\]with a diffusive flux given by
\[F_{diff} = - \rho K_{zz} \frac{\partial q}{\partial z}\]The equation is solved with an implicit scheme assuming that there is no flux at the top and bottom boundaries (evaporation must be treated separately for now).
For the potential temperature mixing, the equation in not completely nrj conserving. One should mix enthalpy as in GCM.
- Parameters:
timestep (
float
) – Time step in seconds.Nlay (
int
) – Number of atmospheric layerst_lay_ov_mu (
array
,np.ndarray
) – Temperatures of the atmospheric layers divided by the molar_mass in kg/molp_lay (
array
,np.ndarray
) – Pressure at the layer centers (Pa)p_lev (
array
,np.ndarray
) – Presure at the Nlay+1 level boundaries (Pa)dmass (
array
,np.ndarray
) – Mass of gas in each layer (kg/m^2)g (
float
) – Gravity (m/s^2)Kzz (
float
) – Eddy mixing coefficient (m^2/s)tracer_array (
array
,np.ndarray (Ntrac
,Nlay)
) – Array containing the mass mixing ratio of all tracers at each layer before the mixing
- Returns:
- new_tracers: array, np.ndarray (Ntrac, Nlay)
Array containing the mass mixing ratio of all tracers at each layer after the mixing
- exo_k.atm_evolution.convection.compute_condensation_numba(timestep, Nlay, tlay, play, cp, Mgas, qarray, idx_vap, idx_cond, thermo_parameters, latent_heating=True, condensation_timestep_reducer=None, verbose=False)[source]
Computes the heating rates needed to bring sursaturated regions back to saturation
- Parameters:
timestep
Nlay (
float
) – Number of layerstlay (
array
,np.ndarray
) – Layer temperaturescp (
float
) – Specific heat capacityMgas (
array
,np.ndarray
) – Layer molar massqarray (
array
,np.ndarray
) – Array of tracer specific concentrationsidx_vap (
int
) – Indices for condensing species in tracer arrayidx_cond (
int
) – Indices for condensing species in tracer arraythermo_parameters (
array
,np.ndarray
) – Array of thermodynamical paramaters for condensing species. See :class:`~exo_k.atm_evolution.condensation.Condensation_Thermodynamical_Parameterslatent_heating (
bool
) – Whether latent heating should be taken into account
- Returns:
- H_cond: array, np.ndarray
Heating rate in W/kg
- exo_k.atm_evolution.convection.moist_convective_adjustment_numba(timestep, Nlay, tlay, play, dmass, cp, Mgas, q_array, i_vap, i_cond, thermo_parameters, moist_inhibition=True, verbose=False)[source]
Computes the heating rates needed to adjust unstable regions of a given atmosphere to a moist adiabat.
Based on Description algorithmique d’un ensemble de paramétrisation physique - phylmd Zhaoxin Laurent LI
- Parameters:
timestep
Nlay (
float
) – Number of layerstlay (
array
,np.ndarray
) – Layer temperaturesplay (
array
) – Pressure at layer centersdmass (
array
,np.ndarray
) – mass of layers in kg/m^2cp (
float
) – specific heat capacity at constant pressureq_array (
array
,np.ndarray
) – mass mixing ratio of tracersi_vap (
int
) – index of vapor tracer in qarrayi_cond (
int
) – index of condensate tracer in qarrayqsat (
array
,np.ndarray
) – Saturation mmr for each layerdqsat_dt (
array
,np.ndarray
) – d qsat / dT in each layerLvap (
array
,np.ndarray
) – Latent heat of vaporization (can have different values in each layer if Lvap=f(T))dlnt_dlnp_moist (
array
,np.ndarray
) – threshold thermal gradient (d ln T / d ln P) for a moist atmosphere computed at the layer centers.q_crit (
array
,np.ndarray
) – Critical mass mixing ratio for the inhibition of moist convection (Eq. 17 of Leconte et al. 2017)
- Returns:
- H_madj: array, np.ndarray
Heating rate in each atmospheric layer (W/kg).
- new_q: array, np.ndarray
tracer mmr array after adjustment.
- new_t: array, np.ndarray
Temperature of layers after adjustment.
- exo_k.atm_evolution.convection.moist_convective_adjustment_cloud_fraction_numba(timestep, Nlay, tlay, play, dmass, cp, Mgas, q_array, i_vap, i_cond, thermo_parameters, moist_inhibition=True, verbose=False, humidity_distribution_width=0.2)[source]
Computes the heating rates needed to adjust unstable regions of a given atmosphere to a moist adiabat.
This approach with a cloud fraction is experimental and should be used with extreme caution. It does not seem to give good results when re-evaporation is very efficient because there is constant re-evaporation and condensation in the convective region.
- Parameters:
timestep
Nlay (
float
) – Number of layerstlay (
array
) – Layer temperaturesplay (
array
) – Pressure at layer centersdmass (
array
) – mass of layers in kg/m^2cp (
float
) – specific heat capacity at constant pressureq_array (
array
) – mass mixing ratio of tracersi_vap (
int
) – index of vapor tracer in qarrayi_cond (
int
) – index of condensate tracer in qarrayqsat (
array
) – Saturation mmr for each layerdqsat_dt (
array
) – d qsat / dT in each layerLvap (
array
) – Latent heat of vaporization (can have different values in each layer if Lvap=f(T))dlnt_dlnp_moist (
array
) – threshold thermal gradient (d ln T / d ln P) for a moist atmosphere computed at the layer centers.q_crit (
array
) – Critical mass mixing ratio for the inhibition of moist convection (Eq. 17 of Leconte et al. 2017)
- Returns:
- H_madj: array
Heating rate in each atmospheric layer (W/kg).
- new_q: array
tracer mmr array after adjustment.
- new_t: array
Temperature of layers after adjustment.
- exo_k.atm_evolution.convection.compute_rainout_numba(timestep, Nlay, tlay, play, dmass, cp, Mgas, qarray, idx_vap, idx_cond, thermo_parameters, evap_coeff, qvap_deep, q_cloud=0.0, latent_heating=False, total_cloud_fraction=1.0, humidity_distribution_width=0.2, verbose=False)[source]
Computes the heating rates needed to adjust unstable regions of a given atmosphere to a moist adiabat.
- Parameters:
timestep
Nlay (
float
) – Number of layerstlay (
array
,np.ndarray
) – Layer temperatures
- exo_k.atm_evolution.convection.molecular_diffusion_numba(timestep, Nlay, p_lay, p_lev, dmass, tlay, mu, g, Dmol, verbose=False)[source]
Solves turbulent diffusion equation:
\[\rho frac{\partialT}{\partial t} = \frac{\partial F_{diff}}{\partial z}\]with a diffusive flux given by
\[F_{diff} = - \rho D_{mol} \frac{\partial T}{\partial z}\]The equation is solved with an implicit scheme assuming that there is no flux at the top and bottom boundaries (evaporation must be treated separately for now).
- Parameters:
timestep (
float
) – Time step in seconds.Nlay (
int
) – Number of atmospheric layerst_lay_ov_mu (
array
,np.ndarray
) – Temperatures of the atmospheric layers divided by the molar_mass in kg/molp_lay (
array
,np.ndarray
) – Pressure at the layer centers (Pa)p_lev (
array
,np.ndarray
) – Presure at the Nlay+1 level boundaries (Pa)dmass (
array
,np.ndarray
) – Mass of gas in each layer (kg/m^2)g (
float
) – Gravity (m/s^2)Dmol (
float
) – molecular diffusion coefficient (m^2/s)
- Returns:
- new_tlay: array, np.ndarray (Nlay)
Array containing the temperature at each layer after the mixing
- exo_k.atm_evolution.convection.DTRIDGL(L, AF, BF, CF, DF)[source]
! GCM2.0 Feb 2003
! DOUBLE PRECISION VERSION OF TRIDGL
DIMENSION AF(L),BF(L),CF(L),DF(L),XK(L) DIMENSION AS(2*L),DS(2*L)
!* THIS SUBROUTINE SOLVES A SYSTEM OF TRIDIAGIONAL MATRIX !* EQUATIONS. THE FORM OF THE EQUATIONS ARE: !* A(I)*X(I-1) + B(I)*X(I) + C(I)*X(I+1) = D(I)
!======================================================================!