Development manual

Contributing

If you are interested in contributing to the model, please contact the developers using the information on the main page.

Coding style

  • Indentation uses 4 spaces, no tabs.
  • Function names should be lowercase, with words separated by underscores .
  • Lines should aim to have a length of no more than 92 characters.
  • All functions should have docstrings, ideally with Arguments and Returns listed.
  • More comments are always better, even if they seem redundant.
  • Use type hinting where possible.
  • Print statements should be made through the logger where possible.
  • The core package code should not contain global variables, except in the phys module.

Code reference

Atmosphere initialisation and variables

Functions from atmosphere.jl.

AGNI.atmosphere.allocate!Method

Allocate atmosphere arrays, prepare spectral files, and final steps.

Will not modify spectral file if stellar_spectrum` is an empty string.

Arguments:

  • atmos::Atmos_t the atmosphere struct instance to be used.
  • stellar_spectrum::String path to stellar spectrum csv file
source
AGNI.atmosphere.calc_layer_props!Method

Calculate properties within each layer of the atmosphere (e.g. density, mmw).

Arguments: - atmos::Atmos_t the atmosphere struct instance to be used. - ignore_errors::Bool do not generate errors from hydrostatic integrator.

Returns: - ok::Bool function result is ok

source
AGNI.atmosphere.calc_observed_rho!Method

Calculate observed radius and bulk density.

This is done at the layer probed in transmission, which is set to a fixed pressure.

Arguments: - atmos::Atmos_t the atmosphere struct instance to be used.

Returns: - transspec_rho::Float64 the bulk density observed in transmission

source
AGNI.atmosphere.calc_profile_cpkc!Method

Calculate specific heat capacity and thermal conductivity for all layers.

Specific heat per unit mass: J K-1 kg-1. Thermal conductivity: W m-1 K-1.

Arguments: - atmos::Atmos_t the atmosphere struct instance to be used.

source
AGNI.atmosphere.calc_profile_density!Method

Calculate the mass-density for all layers.

Requires temperature, pressure, mmw to be already have been set.

Arguments: - atmos::Atmos_t the atmosphere struct instance to be used.

source
AGNI.atmosphere.calc_profile_radius!Method

Calculate radius and gravity for all layers.

Performs hydrostatic integration from the ground upwards. Requires density, temperature, pressure to have already been set.

Arguments: - atmos::Atmos_t the atmosphere struct instance to be used. - ignore_errors::Bool do not generate errors from hydrostatic integrator.

Returns: - ok::Bool function result is ok

source
AGNI.atmosphere.calc_single_cpkc!Method

Calculate specific heat capacity and thermal conductivity of a single layer.

Specific heat per unit mass: J K-1 kg-1. Thermal conductivity: W m-1 K-1.

Arguments: - atmos::Atmos_t the atmosphere struct instance to be used. - idx::Int index of the layer

source
AGNI.atmosphere.calc_single_density!Method

Calculate the mass-density of a single layer.

Requires temperature, pressure, mmw to be already have been set.

Arguments:

  • atmos::Atmos_t the atmosphere struct instance to be used.
  • idx::Int index of the layer
source
AGNI.atmosphere.estimate_photosphere!Method

Estimate photosphere.

Estimates the location of the photosphere by finding the median of the contribution function in each band (0.2 um to 150 um), and then finding the pressure level at which these median values are maximised.

Arguments: - atmos::Atmos_t the atmosphere struct instance to be used.

Returns: - p_ref::Float64 pressure level of photosphere [Pa]

source
AGNI.atmosphere.generate_pgrid!Method

Generate pressure grid.

Almost-equally log-spaced between pboa and pboa. The near-surface layers are smaller than they would be on an equally log-spaced grid, to avoid f numerical weirdness at the bottom boundary.

Arguments:

  • atmos::Atmos_t the atmosphere struct instance to be used.
source
AGNI.atmosphere.integrate_hydrogravMethod

Integrate hydrostatic and gravity equations across a pressure interval.

Uses the classic fourth-order Runge-Kutta method.

Arguments: - r0::Float64 radius at start of interval [m] - g0::Float64 gravity at start of interval [m s-2] - p0::Float64 pressure at start of interval [Pa] - p1::Float64 pressure at end of interval [Pa] - rho::Float64 density throughout interval, constant [kg m-3]

Returns: - r1::Float64 radius at end of interval [m]

source
AGNI.atmosphere.set_tmpl_from_tmp!Method

Set cell-edge temperatures from cell-centre values.

Uses interpolation within the bulk of the column and extrapolation for the topmost edge.

Arguments:

  • atmos::Atmos_t the atmosphere struct instance to be used.
source
AGNI.atmosphere.setup!Method

Set parameters of the atmosphere.

This is used to set up the struct immediately after creation. It must be done before allocate!() is called, and before any RT calculcations are performed.

Arguments:

  • atmos::Atmos_t the atmosphere struct instance to be used.
  • ROOT_DIR::String AGNI root directory.
  • OUT_DIR::String Output directory.
  • spfile::String path to spectral file.
  • instellation::Float64 bolometric solar flux at the top of the atmosphere [W m-2]
  • s0_fact::Float64 scale factor to account for planetary rotation (i.e. S0^*/S0 in Cronin+14)
  • albedo_b::Float64 bond albedo scale factor applied to instellation in order to imitate shortwave reflection
  • zenith_degrees::Float64 angle of radiation from the star, relative to the zenith [degrees].
  • tmp_surf::Float64 effective surface temperature to provide upward longwave flux at the bottom of the atmosphere [K].
  • gravity::Float64 gravitational acceleration at the surface [m s-2].
  • radius::Float64 planet radius at the surface [m].
  • nlev_centre::Int number of model levels.
  • p_surf::Float64 total surface pressure [bar].
  • p_top::Float64 total top of atmosphere pressure [bar].
  • mf_dict::Dict dictionary of VMRs in the format (key,value)=(gas,mf).
  • mf_path::String path to file containing VMRs at each level.

Optional arguments:

  • condensates list of condensates (gas names).
  • surface_material::String surface material (default is "greybody", but can point to file instead).
  • albedo_s::Float64 grey surface albedo used when surface_material="greybody".
  • tmp_floor::Float64 temperature floor [K].
  • C_d::Float64 turbulent heat exchange coefficient [dimensionless].
  • U::Float64 surface wind speed [m s-1].
  • Kzz_floor::Float64 eddy diffusion coefficient, min value [cm2 s-1]
  • tmp_magma::Float64 mantle temperature [K] for sol_type==2.
  • skin_d::Float64 skin thickness [m].
  • skin_k::Float64 skin thermal conductivity [W m-1 K-1].
  • overlap_method::String gaseous overlap scheme (ro: rand overlap, ee: equiv extinct, rorr: ro+resort+rebin).
  • target_olr::Float64 target OLR [W m-2] for sol_type==4.
  • flux_int::Float64 planet's internal flux for sol_type==3.
  • all_channels::Bool use all channels available for RT?
  • flag_rayleigh::Bool include rayleigh scattering?
  • flag_gcontinuum::Bool include generalised continuum absorption?
  • flag_continuum::Bool include continuum absorption?
  • flag_aerosol::Bool include aersols?
  • flag_cloud::Bool include clouds?
  • real_gas::Bool use real gas EOS where possible
  • thermo_functions::Bool use temperature-dependent thermodynamic properties
  • use_all_gases::Bool store information on all supported gases, incl those not provided in cfg

Returns: Nothing

source
AGNI.atmosphere.water_cloud!Method

Manually set water cloud properties at saturated levels.

Uses the default mass mixing ratio, area fraction, and droplet sizes. This function is redundant if condensation is done via chemsitry.handle_saturation!().

Arguments:

  • atmos::Atmos_t the atmosphere struct instance to be used.
source

Energy flux evaluation

Functions from energy.jl.

AGNI.energy.calc_fluxes!Method

Calculate energy flux at each level.

Calculates flux components (radtrans, convection, etc.) and sums them to get total flux. Also updates thermodynamic properties (heat capacity, density, etc.) at each layer.

Arguments:

  • atmos::Atmos_t the atmosphere struct instance to be used.
  • latent::Bool include condensation flux
  • convect::Bool include MLT convection flux
  • sens_heat::Bool include TKE sensible heat transport
  • conduct::Bool include conductive heat transport
  • convect_sf::Float64 scale factor applied to convection fluxes
  • latent_sf::Float64 scale factor applied to phase change fluxes
  • calc_cf::Bool calculate LW contribution function?
  • rainout::Bool allow rainout ( do not reset VMRs to dry values )
source
AGNI.energy.calc_hrates!Method

Calculate heating rates at cell-centres from the total flux.

Requires the total flux to have already been set (atmos.flux_tot). Heating rates are calculated in units of kelvin per day.

Arguments:

  • atmos::Atmos_t the atmosphere struct instance to be used.
source
AGNI.energy.condense_diffuse!Method

Analytical diffusion scheme for condensation and evaporation energy.

Integrates from bottom of model upwards. Based on the amount of phase change at each level, a phase change flux is calculated by assuming a fixed condensation timescale.

Updates fluxes. Requires chemistry.handle_saturation to be called first in the multi-component case.

Arguments:

  • atmos::Atmos_t the atmosphere struct instance to be used.
source
AGNI.energy.convection!Method

Calculate dry convective fluxes using mixing length theory.

Uses the mixing length formulation outlined by Joyce & Tayar (2023), which was also implemented in Lee et al. (2024), and partially outlined in an earlier paper by Robinson & Marley (2014).

https://arxiv.org/abs/2303.09596 https://doi.org/10.1093/mnras/stae537 https://ui.adsabs.harvard.edu/abs/1962JGR....67.3095B/abstract

Convective energy transport fluxes are calculated at every level edge, just like the radiative fluxes. This is not compatible with moist convection. By using MLT to parameterise convection, we can also calculate Kzz directly.

The mixing length is set to asymptotically approach H (for z>>H) or z (for z<H) as per Blackadar (1962). Alternatively, it can be set equal to H.

The scale height is formulated as: Hp = P / (ρ g) The adiabatic lapse rate is formulated as: ∇_ad = dln(T)/dln(P) = (P/T)*(dT/dP) = (P/T)*(1/[ρ c_p]) for an ideal gas, this becomes: ∇_ad = R / (μ c_p)

Arguments:

  • atmos::Atmos_t the atmosphere struct instance to be used.
  • pmin::Float64 pressure [bar] below which convection is disabled
  • mltype::Int mixing length value (1: scale height, 2: asymptotic)
source
AGNI.energy.fill_Kzz!Method

Fill Kzz values for the entire profile.

This function is called after the convection scheme has been run.

The Kzz value in the convective region (and below) are set equal to the maximum value in the convective region, as calculated by MLT. The value increases with power-law scaling with pressure in the stratosphere.

Arguments:

  • atmos::Atmos_t the atmosphere struct instance to be used.
source
AGNI.energy.radtrans!Method

Calculate radiative fluxes using SOCRATES.

Uses the configuration inside the atmos struct. Can either do LW or SW calculation as required. Imports SOCRATES wrapper from the atmosphere module, rather than loading it twice.

Arguments:

  • atmos::Atmos_t the atmosphere struct instance to be used.
  • lw::Bool longwave calculation? Else: shortwave
  • calc_cf::Bool=false also calculate contribution function?
source

Numerical solver

Functions from solver.jl.

AGNI.solver.gs_searchMethod

Golden section search algorithm

Minimises a function f between the bounds a and b. The function f must return a positive scalar.

Arguments:

  • f Function to be minimised
  • a Initial bracket, lower value
  • b Initial bracket, upper valuie
  • dxtol Convergence: exit when bracket is smaller than this size
  • atol Convergence: absolute tolerance on minimum
  • rtol Convergence: relative tolerance on minimum
  • max_steps Maximum number of iterations
  • warnings Print warning if search does not converge before max_steps` are taken.
source
AGNI.solver.solve_energy!Method

Obtain radiative-convective equilibrium using a matrix method.

Solves the non-linear system of equations defined by the flux field divergence, minimising flux loss across a cell by iterating the temperature profile.

Arguments:

  • atmos::Atmos_t the atmosphere struct instance to be used.
  • sol_type::Int solution type, 1: tmpsurf | 2: skin | 3: fluxint | 4: tgt_olr
  • chem_type::Int chemistry type (see wiki)
  • convect::Bool include convection
  • sens_heat::Bool include sensible heating at the surface
  • conduct::Bool include conductive heat transport within the atmosphere
  • latent::Bool include latent heat exchange (condensation/evaporation)
  • rainout::Bool allow rainout (phase change impacts mixing ratios, not just energy fluxes)
  • dx_max::Float64 maximum step size [K]
  • max_steps::Int maximum number of solver steps
  • max_runtime::Float64 maximum runtime in wall-clock seconds
  • fdw::Float64 finite difference: relative width (dx/x) of the "difference"
  • fdc::Bool finite difference: ALWAYS use central difference?
  • fdo::Int finite difference: scheme order (2nd or 4th)
  • method::Int numerical method (1: Newton-Raphson, 2: Gauss-Newton, 3: Levenberg-Marquardt)
  • ls_method::Int linesearch algorithm (0: None, 1: golden, 2: backtracking)
  • easy_start::Bool improve convergence by introducing convection and phase change gradually
  • perturb_all::Bool always recalculate entire Jacobian matrix? Otherwise updates columns only as required
  • ls_increase::Bool factor by which the cost can increase from last step before triggering linesearch
  • detect_plateau::Bool assist solver when it is stuck in a region of small dF/dT
  • modplot::Int iteration frequency at which to make plots
  • save_frames::Bool save plotting frames
  • modprint::Int iteration frequency at which to print info
  • plot_jacobian::Bool plot jacobian too?
  • conv_atol::Float64 convergence: absolute tolerance on per-level flux deviation [W m-2]
  • conv_rtol::Float64 convergence: relative tolerance on per-level flux deviation [dimensionless]

Returns: Nothing

source

Plotting functions and utilities

Functions from plotting.jl.

AGNI.plotting.plot_contfunc2Method

Plot normalised contribution function (per band)

The data displayed in this plot are fine, but the x-axis ticks are labelled incorrectly by the plotting library. I don't know why this is.

source

File I/O modules

Functions from save.jl and load.jl.

AGNI.save.write_fluxesMethod

Write cell-edge energy fluxes to a CSV file

Arguments

  • atmos::Atmos_t Atmosphere object
  • fname::String Filename to write to
source
AGNI.save.write_ncdfMethod

Write verbose atmosphere data to a NetCDF file

Arguments

  • atmos::Atmos_t Atmosphere object
  • fname::String Filename to write to
source
AGNI.save.write_profileMethod

Write {Pressure, Temperature, Radius} profile to a CSV file

Arguments

  • atmos::Atmos_t Atmosphere object
  • fname::String Filename to write to
source
AGNI.load.load_ncdf!Method

Load atmosphere data from a NetCDF file.

This must be called after setup! and allocate! have been called on this atmos struct. Does not overwrite energy fluxes within the atmospherem, or other wavelength-dependent properties. Only composition, temperatures, layer properties, and boundary conditions are updated by this function.

Arguments:

  • atmos::Atmos_t the atmosphere struct instance to be used.
  • fname::String path to the NetCDF file

Returns:

  • ok::Bool loaded ok?
source