exojax.atm package

Submodules

exojax.atm.amclouds module

Ackerman and Marley 2001 cloud model.

  • based on Ackerman and Marley (2001) ApJ 556, 872, hereafter AM01

exojax.atm.amclouds.compute_cloud_base_pressure(pressure, saturation_pressure, vmr_vapor)

computes cloud base pressure from an intersection of a T-P profile and Psat(T) curves :param pressure: pressure array :param saturation_presure: saturation pressure arrau :param vmr_vapor: volume mixing ratio (VMR) for vapor

Returns:

cloud base pressure

exojax.atm.amclouds.compute_cloud_base_pressure_index(pressure, saturation_pressure, vmr_vapor)

computes cloud base pressure from an intersection of a T-P profile and Psat(T) curves :param pressure: pressure array :param saturation_presure: saturation pressure arrau :param vmr_vapor: volume mixing ratio (VMR) for vapor

Returns:

cloud base pressure index

Return type:

int

exojax.atm.amclouds.find_rw(rarr, terminal_velocity, Kzz_over_L)

finding rw from rarr and terminal velocity array.

Parameters:
  • rarr – particle radius array (cm)

  • terminal_velocity – terminal velocity (cm/s)

  • Kzz_over_L – Kzz/L in Ackerman and Marley 2001

Returns:

rw in Ackerman and Marley 2001

exojax.atm.amclouds.geometric_radius(rg, sigmag)

computes the paritculate geometric radius

Parameters:
  • rg (float) – rg parameter in lognormal distribution in cgs

  • sigmag (float) – sigma_g parameter in lognormal distribution

Note

The cross section is given by $S = Q_e pi r_geo^2$

Returns:

geometric radius in cgs

Return type:

_type_

exojax.atm.amclouds.get_rg(rw, fsed, alpha, sigmag)

compute rg of the lognormal size distribution defined by (9) in AM01. The computation is based on (13) in AM01.

Parameters:
  • rw – rw (cm)

  • fsed – fsed

  • alpha – power of the condensate size distribution

  • sigmag – sigmag parameter (geometric standard deviation) in the lognormal distribution of condensate size, defined by (9) in AM01, must be sigmag > 1

Returns

rg: rg parameter in the lognormal distribution of condensate size, defined by (9) in AM01

exojax.atm.amclouds.get_rw(terminal_velocity, Kzz, L, rarr)

compute rw in AM01 implicitly defined by (11)

Parameters:
  • vfs – terminal velocity (cm/s)

  • Kzz – diffusion coefficient (cm2/s)

  • L – typical convection scale (cm)

  • rarr – condensate scale array

Returns:

rw (cm) in AM01. i.e. condensate size that balances an upward transport and sedimentation

Return type:

rw

exojax.atm.amclouds.mixing_ratio_cloud_pressure(pressure, cloud_base_pressure, fsed, mr_cloud_base, kc)

mol mixing ratio of clouds based on AM01 a given single pressure.

Parameters:
  • pressure (float) – pressure (bar) where we want to compute VMR of clouds

  • cloud_base_pressure – cloud base pressure (bar)

  • fsed – fsed

  • mr_cloud_base – mass mixing ratio (MMR) or mol mixing ratio of condensate at cloud base

  • kc – constant ratio of condenstates to total mixing ratio

Returns:

mol mixing ratio of condensates

exojax.atm.amclouds.mixing_ratio_cloud_profile(pressures, cloud_base_pressure, fsed, mr_cloud_base, kc=1)

volume mixing ratio of clouds based on AM01 given pressure.

Parameters:
  • pressures – pressure array (Nlayer) (bar) where we want to compute VMR of clouds

  • cloud_base_pressure – cloud base pressure (bar)

  • fsed – fsed

  • mr_cloud_base – mass mixing ratio (MMR) or mol mixing ratio of condensate at cloud base

  • kc – constant ratio of condenstates to total mixing ratio

Returns:

VMR of condensates

exojax.atm.atmconvert module

exojax.atm.atmconvert.mmr_to_vmr(mmr, molecular_mass, mean_molecular_weight)

convert mass mixing ratio (mmr) to volume mixing ratio (vmr)

Parameters:
  • mmr (float or array) – mass mixing ratio(s)

  • molecular_mass (float or array) – molecular mass(es)

  • mean_molecular_weight (float) – mean molecular weight

Returns:

volume mixing ratio(s)

Return type:

float

exojax.atm.atmconvert.vmr_to_mmr(vmr, molecular_mass, mean_molecular_weight)

convert volume mixing ratio (vmr) to mass mixing ratio (mmr)

Parameters:
  • vmr (float or array) – volume mixing ratio(s)

  • molecular_mass (float or array) – molecular mass(es)

  • mean_molecular_weight (float) – mean molecular weight

Returns:

mass mixing ratio(s)

Return type:

float

exojax.atm.atmphys module

Atmospheric MicroPhysics (amp) class

exojax.atm.atmprof module

Atmospheric profile function.

exojax.atm.atmprof.Teff2Tirr(Teff, Tint)

Tirr from effective temperature and intrinsic temperature.

Parameters:
  • Teff – effective temperature

  • Tint – intrinsic temperature

Returns:

iradiation temperature

Return type:

Tirr

Note

Here we assume A=0 (albedo) and beta=1 (fully-energy distributed)

exojax.atm.atmprof.Teq2Tirr(Teq)

Tirr from equilibrium temperature and intrinsic temperature.

Parameters:

Teq – equilibrium temperature

Returns:

iradiation temperature

Return type:

Tirr

Note

Here we assume A=0 (albedo) and beta=1 (fully-energy distributed)

exojax.atm.atmprof.atmprof_Guillot(pressures, gravity, kappa, gamma, Tint, Tirr, f=0.25)

Notes

Guillot (2010) Equation (29)

Parameters:
  • pressures – pressure array (bar)

  • gravity – gravity (cm/s2)

  • kappa – thermal/IR opacity (kappa_th in Guillot 2010)

  • gamma – ratio of optical and IR opacity (kappa_v/kappa_th), gamma > 1 means thermal inversion

  • Tint – temperature equivalence of the intrinsic energy flow

  • Tirr – temperature equivalence of the irradiation

  • point (f = 1 at the substellar) – and f = 1/4 for an averaging over the whole planetary surface

  • average (f = 1/2 for a day-side) – and f = 1/4 for an averaging over the whole planetary surface

Returns:

temperature profile

Return type:

array

exojax.atm.atmprof.atmprof_gray(pressures, gravity, kappa, Tint)
Parameters:
  • pressures (1D array) – pressure array (bar)

  • gravity (float) – gravity (cm/s2)

  • kappa – infrared opacity

  • Tint – temperature equivalence of the intrinsic energy flow

Returns:

temperature profile

Return type:

array

exojax.atm.atmprof.atmprof_powerlow(pressures, T0, alpha)

powerlaw temperature profile

Parameters:
  • pressures – pressure array (bar)

  • T0 (float) – T at P=1 bar in K

  • alpha (float) – powerlaw index

Returns:

temperature profile

Return type:

array

exojax.atm.atmprof.gh_product(temperature, mean_molecular_weight)

product of gravity and pressure scale height

Parameters:
  • temperature – isothermal temperature (K)

  • mean_molecular_weight – mean molecular weight

Returns:

gravity x pressure scale height cm2/s2

exojax.atm.atmprof.normalized_layer_height(temperature, pressure_decrease_rate, mean_molecular_weight, radius_btm, gravity_btm)

compute normalized height/radius at the upper boundary of the atmospheric layer, neglecting atmospheric mass, examining non-constant gravity.

Note

This method computes the height of the atmospheric layers taking the effect of the decrease of gravity (i.e. $ propto 1/r^2 $) into account.

Parameters:
  • temperature (1D array) – temperature profile (K) of the layer, (Nlayer, from atmospheric top to bottom)

  • pressure_decrease_rate – pressure decrease rate of the layer (k-factor; k < 1) pressure[i-1] = pressure_decrease_rate*pressure[i]

  • mean_molecular_weight (1D array) – mean molecular weight profile, (Nlayer, from atmospheric top to bottom)

  • radius_btm (float) – radius (cm) at the lower boundary of the bottom layer, R0 or r_N

  • gravity_btm (float) – gravity (cm/s2) at the lower boundary of the bottom layer, g_N

Returns:

layer height normalized by radius_btm starting from top atmosphere 1D array (Nlayer) : radius at lower bondary normalized by radius_btm starting from top atmosphere

Return type:

1D array (Nlayer)

exojax.atm.atmprof.pressure_boundary_logspace(pressures, pressure_decrease_rate, reference_point=0.5, numpy=False)

computes pressure at the boundary of the layers (Nlayer + 1)

Parameters:
  • pressures (_type_) – representative pressure (output of pressure_layer_logspace)

  • pressure_decrease_rate – pressure decrease rate of the layer (k-factor; k < 1) pressure[i-1] = pressure_decrease_rate*pressure[i]

  • reference_point (float) – reference point in a layer (0-1). Center:0.5, lower boundary:1.0, upper boundary:0

  • numpy – if True use numpy array instead of jnp array

Returns:

pressure at the boundary (Nlayer + 1)

Return type:

_type_

exojax.atm.atmprof.pressure_layer_logspace(log_pressure_top=- 8.0, log_pressure_btm=2.0, nlayer=20, mode='ascending', reference_point=0.5, numpy=False)

Pressure layer evenly spaced in logspace, i.e. logP interval is constant

Parameters:
  • log_pressure_top – log10(P[bar]) at the top layer

  • log_pressure_btm – log10(P[bar]) at the bottom layer

  • nlayer – the number of the layers

  • mode – ascending or descending

  • reference_point (float) – reference point in a layer (0-1). Center:0.5, lower boundary:1.0, upper boundary:0

  • numpy – if True use numpy array instead of jnp array

Returns:

representative pressures (array) of the layers delta_pressures: delta pressure layer, the old name is dParr pressure_decrease_rate: pressure decrease rate of the layer (k-factor; k < 1) pressure[i-1] = pressure_decrease_rate*pressure[i]

Return type:

pressures

Note

d logP is constant using this function.

exojax.atm.atmprof.pressure_lower_logspace(pressures, pressure_decrease_rate, reference_point=0.5)

computes pressure at the lower point of the layers

Parameters:
  • pressures (_type_) – representative pressure (output of pressure_layer_logspace)

  • pressure_decrease_rate – pressure decrease rate of the layer (k-factor; k < 1) pressure[i-1] = pressure_decrease_rate*pressure[i]

  • reference_point (float) – reference point in a layer (0-1). Center:0.5, lower boundary:1.0, upper boundary:0

Returns:

pressure at the lower point (underline{P}_i)

Return type:

_type_

exojax.atm.atmprof.pressure_scale_height(gravity, T, mean_molecular_weight)

pressure scale height assuming an isothermal atmosphere.

Parameters:
  • gravity – gravity acceleration (cm/s2)

  • T – isothermal temperature (K)

  • mean_molecular_weight – mean molecular weight

Returns:

pressure scale height (cm)

exojax.atm.atmprof.pressure_upper_logspace(pressures, pressure_decrease_rate, reference_point=0.5)

computes pressure at the upper point of the layers

Parameters:
  • pressures (_type_) – representative pressure (output of pressure_layer_logspace)

  • pressure_decrease_rate – pressure decrease rate of the layer (k-factor; k < 1) pressure[i-1] = pressure_decrease_rate*pressure[i]

  • reference_point (float) – reference point in a layer (0-1). Center:0.5, lower boundary:1.0, upper boundary:0

Returns:

pressure at the upper point (overline{P}_i)

Return type:

_type_

exojax.atm.condensate module

exojax.atm.condensate.condensate_density_liquid_ammonia(T)

condensate density of liquid ammonia as a function of T from Lange’s Handbook of Chemistry, 10th ed. page 1451 and 1468

Parameters:

T (_type_) – temperature in Kelvin

Returns:

liquid ammonia density in g/cm3

exojax.atm.condensate.read_liquid_ammonia_density()

read liquid ammonia density file

Returns:

temperature in K liquid ammonia density in g/cm3

exojax.atm.idealgas module

Functions about ideal gas.

exojax.atm.idealgas.number_density(Parr, Tarr)

number density of ideal gas in cgs.

Parameters:
  • Parr – pressure array (bar)

  • Tarr – temperature array (K)

Returns:

number density (1/cm3)

exojax.atm.lorentz_lorenz module

exojax.atm.lorentz_lorenz.refractive_index_Lorentz_Lorenz(polarizability, number_density)

Refractive index using Lorentz-Lorenz relataion

Notes

See D.25 in Liou 2002, for instance.

Parameters:
  • polarizability – polarizability (cm3)

  • number_density – number density of molecule (cm-3)

Returns:

refractive index

exojax.atm.mixratio module

Mixing ratio definition

exojax.atm.mixratio.mmr2vmr(mmr, molecular_mass, mean_molecular_weight)

MMR (Mass Mixing Ratio) to VMR (Volume Mixing Ratio or Mol Mixing ratio)

Parameters:
  • mmr (_type_) – Mass Mixing Ratio

  • molecular_mass (_type_) – Molecular Mass (i.e. utils.molinfo.molmass_isotope(“H2O”) for instance)

  • mean_molecular_weight (_type_) – Mean Molecular Weight

Returns:

Volume Mixing Ratio or Mol Mixing Ratio

Return type:

_type_

exojax.atm.mixratio.vmr2mmr(vmr, molecular_mass, mean_molecular_weight)

VMR (Volume Mixing Ratio or Mol Mixing ratio) to MMR (Mass Mixing Ratio)

Parameters:
  • vmr (_type_) – Volume Mixing Ratio

  • molecular_mass (_type_) – Molecular Mass (i.e. utils.molinfo.molmass_isotope(“H2O”) for instance)

  • mean_molecular_weight (_type_) – Mean Molecular Weight

Returns:

Mass Mixing Ratio

Return type:

_type_

exojax.atm.polarizability module

Gas Polarizability

Notes

Originally taken from PICASO/GPLv3 picaso/rayleigh.py polarizabilities are mainly taken from CRC handbook of chemistry and physics vol. 95 unit=cm3 H3+ taken from Kawaoka & Borkman, 1971 Number density at reference conditions of refractive index measurements i.e. number density of the ideal gas at T=273.15K (=0 C) and P=1atm [cm-2], as Patm*bar_cgs / (kB * 273.15) http://refractiveindex.info n_ref_refractive = 2.6867810458916872e+19

exojax.atm.psat module

Saturation Vapor Pressure.

exojax.atm.psat.psat_Fe_AM01(T)

Saturation Vapor Pressure for iron

Note

Taken from Ackerman and Marley 2001 Appendix A (A3), originally from Barshay and Lewis (1976)

Parameters:

T – temperature (K)

Returns:

saturation vapor pressure (bar)

exojax.atm.psat.psat_ammonia_AM01(T)

Saturation Vapor Pressure for Ammonia

Note

Taken from Ackerman and Marley 2001 Appendix A (A4) see also their errata.

Parameters:

T – temperature (K)

Returns:

saturation vapor pressure (bar)

exojax.atm.psat.psat_enstatite_AM01(T)

Saturation Vapor Pressure for Enstatite (MgSiO3)

Note

Taken from Ackerman and Marley 2001 Appendix A (A4) see also their errata, originally from Barshay and Lewis (1976)

Parameters:

T – temperature (K)

Returns:

saturation vapor pressure (bar)

exojax.atm.psat.psat_water_AM01(T)

Saturation Vapor Pressure for water (but for T<1048K)

Note

Taken from Ackerman and Marley 2001 Appendix A (A4) see also their errata and updated Buck 1981, 1996

Parameters:

T – temperature (K)

Returns:

saturation vapor pressure (bar)

exojax.atm.psat.psat_water_Magnus(T)

Saturation Vapor Pressure for water (Magnus, or August-Roche-Magnus)

Parameters:

T – temperature (K)

Returns:

saturation vapor pressure (bar)

exojax.atm.simple_clouds module

cloud opacity.

exojax.atm.simple_clouds.powerlaw_clouds(nus, kappac0=0.01, nuc0=28571.0, alphac=1.0)

power-law cloud model.

Parameters:
  • kappac0 – opacity (cm2/g) at nuc0

  • nuc0 – wavenumber for kappac0

  • alphac – power

Returns:

cross section (cm2)

Note

alphac = - gamma of the definition in petitRadtrans. Also, default nuc0 corresponds to lambda0=0.35 um.

exojax.atm.viscosity module

Viscosity of droplets.

exojax.atm.viscosity.calc_vfactor(atm='H2', LJPparam=None)
Parameters:
  • atm – molecule consisting of atmosphere, “H2”, “O2”, and “N2”

  • LJPparam – Custom Lennard-Jones Potential Parameters (d (cm) and epsilon/kB)

Returns:

dynamic viscosity factor for Rosner eta = viscosity*T**0.66 applicable tempature range (K,K)

Return type:

vfactor

Note

The dynamic viscosity is from the Rosner book (3-2-12) and caption in p106 Hirschfelder et al. (1954) within Trange.

exojax.atm.viscosity.eta_Rosner(T, vfactor)

dynamic viscocity by Rosner (2000)

Parameters:
  • T – temperature (K)

  • vfactor – vfactor

Returns:

dynamic viscosity (g/s/cm)

exojax.atm.viscosity.eta_Rosner_H2(T)

dynamic viscocity of the H2 atmosphere by Rosner (2000)

Parameters:

T – temperature (K) (applicable from 179 to 11940K)

Returns:

dynamic viscosity (g/s/cm)

exojax.atm.viscosity.get_LJPparam()

Lennard-Jones Potential Parameters.

Returns:

Dict for Lennard-Jones Potential Parameters (d (cm)) LJPparam_epsilon_per_kB: Dict for Lennard-Jones Potential Parameters (epsilon/kB)

Return type:

LJPparam_d

Note

Lennard-Jones Potential Parameters (LJPparam) were taken from p107, Table 3.2.1 of Transport process in chemically reacting flow systems by Daniel E. Rosner, originally from Svehla (1962).

exojax.atm.vterm module

Terminal velocity of cloud particles.

Note

The code in this module is based on Hans R Pruppacher and James D Klett. Microstructure of atmospheric clouds and precipitation and Akerman and Marley 2001.

exojax.atm.vterm.Ndavies(r, g, eta, drho, rho)

Davies (Best) number.

Parameters:
  • r – particle size (cm)

  • g – gravity (cm/s2)

  • eta – dynamic viscosity (g/s/cm)

  • drho – density difference between condensates and atmosphere (g/cm3)

  • rho – atmosphere density (g/cm3)

Returns:

Davies number (Best Number)

exojax.atm.vterm.terminal_velocity(r, gravity, dynamic_viscosity, rho_cloud, rho_atm, Nkn=0.0)

computes terminal velocity in a wide particles size range.

Parameters:
  • r – particle size (cm)

  • g – gravity (cm/s2)

  • dynamic_viscosity – dynamic viscosity (g/s/cm)

  • rho_cloud – condensate density (g/cm3)

  • rho_atm – atmosphere density (g/cm3)

  • Nkn – Knudsen number

Returns:

terminal velocity (cm/s)

Example

>>> #terminal velocity at T=300K, for Earth atmosphere/gravity.
>>> g=980.
>>> drho=1.0
>>> rho=1.29*1.e-3 #g/cm3
>>> vfactor,Tr=vc.calc_vfactor(atm="Air")
>>> eta=vc.eta_Rosner(300.0,vfactor)
>>> r=jnp.logspace(-5,0,70)
>>> terminal_velocity(r,g,eta,drho,rho) #terminal velocity (cm/s)
exojax.atm.vterm.vf_largeNre(r, g, eta, drho, rho)

terminal velocity ( Reynolds number > 500, Davies number >10**5 )

Parameters:
  • r – particle size (cm)

  • g – gravity (cm/s2)

  • eta – dynamic viscosity (g/s/cm)

  • drho – density difference between condensates and atmosphere (g/cm3)

  • rho – atmosphere density (g/cm3)

Returns:

terminal velocity (cm/s)

exojax.atm.vterm.vf_midNre(r, g, eta, drho, rho)

terminal velocity (2 < Reynolds number < 500, 42 < Davies number < 10**5)

Parameters:
  • r – particle size (cm)

  • g – gravity (cm/s2)

  • eta – dynamic viscosity (g/s/cm)

  • drho – density difference between condensates and atmosphere (g/cm3)

  • rho – atmosphere density (g/cm3)

Returns:

terminal velocity (cm/s)

exojax.atm.vterm.vf_stokes(r, g, eta, drho, Nkn=0.0)

terminal velocity of Stokes flow (Reynolds number < 2, Davies number < 42)

Parameters:
  • r – particle size (cm)

  • g – gravity (cm/s2)

  • eta – dynamic viscosity (g/s/cm)

  • drho – density difference between condensates and atmosphere (g/cm3)

  • Nkn – Knudsen number

Returns:

terminal velocity (cm/s)

Note

(1.0+1.255*Nkn) is the Cunningham factor

Note

See also (10-138) p415 in Hans R Pruppacher and James D Klett. Microstructure of atmospheric clouds and precipitation. In Microphysics of clouds and precipitation, pages 10–73. Springer, 2010. Equation (B1) in Appendix B of Ackerman and Marley 2001.

Module contents