Software manual
Contributing
If you are interested in contributing to the model, please see the contributing page
Code reference (API)
Atmosphere initialisation and variables
Functions from atmosphere.jl.
AGNI.atmosphere — Module
Main atmosphere module for storing model data
This module primarily contains the Atmos_t struct which does most of the heavy lifting in AGNI. There are also functions for setting up and configuring the struct.
Also includes hydrostatic integrator and ways to estimate observable quantities.
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. Will treat star as blackbody with photospheric effective temperature of stellar_Teff if the parameter stellar_spectrum has value of "blackbody".
Arguments:
atmos::Atmos_tthe atmosphere struct instance to be used.stellar_spectrum::Stringpath to stellar spectrum csv filestellar_Teff::Float64star effective temperature if blackbodycheck_safe_gas::Boolrequire that there be at least one 'safe' gas in the mix
AGNI.atmosphere.calc_layer_props! — Method
Calculate properties within each layer of the atmosphere (e.g. density, mmw).
Function will return false if hydrostatic calculcation fails. This is usually when the atmosphere becomes unbound.
Arguments: - atmos::Atmos_t the atmosphere struct instance to be used.
Returns: - ok::Bool function result is ok
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
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.
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.
AGNI.atmosphere.calc_profile_mmw! — Method
Calculate mean molecular weight for all layers.
MMW stored as kg/mol.
Arguments: - atmos::Atmos_t the atmosphere struct instance to be used.
AGNI.atmosphere.calc_profile_radius! — Method
Calculate radii, gravities, and masses for all layers.
Performs hydrostatic integration from the ground upwards. Requires density, temperature, pressure to have already been set.
Does not account for surface ocean height.
Arguments:
atmos::Atmos_tthe atmosphere struct instance to be used.
Returns:
bound::Boolatmosphere is strongly bound by gravity
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
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_tthe atmosphere struct instance to be used.idx::Intindex of the layer
AGNI.atmosphere.estimate_Ra! — Method
Estimate a diagnostic Rayleigh number in each layer.
Assuming that the Rayleigh number scales like Ra ~ (wλ/κ)^(1/β) Where κ is the thermal diffusivity and β is the convective beta parameter.
This quantity must be taken lightly.
Arguments:
atmos::Atmos_tthe atmosphere struct instance to be used.
AGNI.atmosphere.estimate_convective_zone — Method
Get pressure at top and bottom of convective zone
Arguments: - atmos::Atmos_t the atmosphere struct instance to be used.
Returns: - p_top::Float64 pressure [Pa] at top of convective zone - p_bot::Float64 pressure [Pa] at bottom of convective zone
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]
AGNI.atmosphere.estimate_timescale_conv! — Method
Estimate a diagnostic convective timescale in each layer.
This quantity must be taken lightly.
Arguments:
atmos::Atmos_tthe atmosphere struct instance to be used.
AGNI.atmosphere.estimate_timescale_rad! — Method
Estimate a diagnostic radiative timescale in each layer.
This quantity must be taken lightly.
Arguments:
atmos::Atmos_tthe atmosphere struct instance to be used.
AGNI.atmosphere.generate_pgrid! — Method
Generate pressure grid.
Almost equally log-spaced between pboa and pboa. The near-boundary layers are smaller than they would be on an equally log-spaced grid, to avoid numerics.
Arguments:
atmos::Atmos_tthe atmosphere struct instance to be used.
AGNI.atmosphere.get_interleaved_ptr — Method
Get interleaved cell-centre and cell-edge PTR grid.
Arguments:
atmos::Atmos_tthe atmosphere struct instance to be used.
Returns:
arr_P::Arraypressure array [Pa]arr_T::Arraytemperature array [K]arr_R::Arrayradius array [m]
AGNI.atmosphere.integ_hydrograv — Method
Integrate hydrostatic and gravity equations across a pressure interval.
Uses the classic fourth-order Runge-Kutta method.
Arguments:
r0::Float64radius at start of interval [m]g0::Float64gravity at start of interval [m s-2]m0::Float64mass enc at start of interval [kg]p0::Float64pressure at start of interval [Pa]p1::Float64pressure at end of interval [Pa]rho::Float64density throughout interval, constant [kg m-3]n::Intnumber of steps for integration (n >= 2)
Returns:
rj::Float64radius at end of interval [m]gj::Float64gravity at end of interval [kg]mj::Float64mass enc at end of interval [kg]
AGNI.atmosphere.make_transparent! — Method
Set atmosphere properties such that it is effectively transparent.
This will modify the surface pressure and disable gas opacity in SOCRATES. These changes cannot be reversed directly. To undo them, it is best to create and allocate a new atmosphere struct.
Arguments:
atmos::Atmos_tthe atmosphere struct instance to be used.
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_tthe atmosphere struct instance to be used.
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_tthe atmosphere struct instance to be used.ROOT_DIR::StringAGNI root directory.OUT_DIR::StringOutput directory.spfile::Stringpath to spectral file ('grey' => use grey gas)instellation::Float64bolometric solar flux at the top of the atmosphere [W m-2]s0_fact::Float64scale factor to account for planetary rotation (i.e. S0^*/S0 in Cronin+14)albedo_b::Float64bond albedo scale factor applied to instellation in order to imitate shortwave reflectionzenith_degrees::Float64angle of radiation from the star, relative to the zenith [degrees].tmp_surf::Float64effective surface temperature to provide upward longwave flux at the bottom of the atmosphere [K].gravity::Float64gravitational acceleration at the surface [m s-2].radius::Float64planet radius at the surface [m].nlev_centre::Intnumber of model levels.p_surf::Float64total surface pressure [bar].p_top::Float64total top of atmosphere pressure [bar].mf_dict::Dictdictionary of VMRs in the format (key,value)=(gas,mf).mf_path::Stringpath to file containing VMRs at each level.
Optional arguments:
IO_DIR::Stringdirectory used for fast file operations.condensateslist of condensates (gas names).metallicities::Dictdictionary of elemental metallicities (mass ratio rel to hydrogen)surface_material::Stringsurface material (default is "greybody", but can point to file instead).albedo_s::Float64grey surface albedo used whensurface_material="greybody".tmp_floor::Float64temperature floor [K].surf_roughness::Float64surface roughness length scale [m]surf_windspeed::Float64surface wind speed [m s-1].Kzz_floor::Float64min eddy diffusion coefficient, cgs units [cm2 s-1]mlt_asymptotic::Boolmixing length scales asymptotically, but ~0 near groundmlt_criterion::CharMLT stability criterion. Options: (s)chwarzschild, (l)edoux.tmp_magma::Float64mantle temperature [K] for sol_type==2.skin_d::Float64skin thickness [m].skin_k::Float64skin thermal conductivity [W m-1 K-1].overlap_method::Stringgaseous overlap scheme (ro: rand overlap, ee: equiv extinct, rorr: ro+resort+rebin).target_olr::Float64target OLR [W m-2] for sol_type==4.flux_int::Float64planet's internal flux for sol_type==3.all_channels::Booluse all channels available for RT?flag_rayleigh::Boolinclude rayleigh scattering?flag_gcontinuum::Boolinclude generalised continuum absorption?flag_continuum::Boolinclude continuum absorption?flag_aerosol::Boolinclude aersols?flag_cloud::Boolinclude clouds?phs_timescale::Float64phase change timescale [s]evap_efficiency::Float64re-evaporatione efficiency compared to saturating amountreal_gas::Booluse real gas EOS where possiblethermo_functions::Booluse temperature-dependent thermodynamic propertiesuse_all_gases::Boolstore information on all supported gases, incl those not provided in cfgcheck_integrity::Boolconfirm integrity of thermo files using their checksumκ_grey_lw::Float64gas opacity when using grey-gas RT scheme, longwaveκ_grey_sw::Float64gas opacity when using grey-gas RT scheme, shortwavefastchem_work::Stringworking directory for fastchemfastchem_floor::Float64temperature floor on profile provided to fastchemfastchem_maxiter_chem::Float64maximum chemical iterations allowed by fastchemfastchem_maxiter_solv:Float64maximum solver iterations allowed by fastchemfastchem_xtol_chem::Float64solution tolerance required of fastchem (chemical)fastchem_xtol_elem::Float64solution tolerance required of fastchem (elemental)rfm_parfile::Stringpath to HITRAN-formatted .par file provided to RFM
Returns: Nothing
Energy flux calculations
Functions from energy.jl.
AGNI.energy._radtrans_greygas! — Method
Solve RT using double grey-gas formulation
Simple two-stream double grey RT solver which integrates fluxes from the TOA and BOA.
Uses two opacity values to represent the LW and SW components of the flux field.
Loosely following this tutorial, which is based on Pierrehumbert (2010). https://brian-rose.github.io/ClimateLaboratoryBook/courseware/radiative-transfer/
Arguments:
atmos::Atmos_tthe atmosphere struct instance to be used.
AGNI.energy._radtrans_socrates! — Method
Solve radiative transfer using SOCRATES
Imports SOCRATES wrapper from the atmosphere module, rather than loading it twice.
Arguments:
atmos::Atmos_tthe atmosphere struct instance to be used.lw::BoolTrue: longwave calculation. False: shortwave calculation.
Optional arguments:
calc_cf::Boolalso calculate contribution function?gauss_ir::Boolusing gaussian angular integration in IR, otherwise uses two-stream approximationrescale_pf::Boolperform rescaling on phase function
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.
Assumes that chemistry functions have already been called, if wanted. Does not call fastchem here.
Arguments:
atmos::Atmos_tthe atmosphere struct instance to be used.
Optional arguments:
radiative::Boolinclude radiation fluxeslatent_heat::Boolinclude condensation fluxconvective::Boolinclude MLT convection fluxsens_heat::Boolinclude TKE sensible heat fluxconductive::Boolinclude conductive heat fluxadvective::Boolinclude advective heat fluxconvect_sf::Float64scale factor applied to convection fluxeslatent_sf::Float64scale factor applied to phase change fluxescalc_cf::Boolcalculate LW contribution function?
Returns:
ok::Boolcalculation performed ok?
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_tthe atmosphere struct instance to be used.
AGNI.energy.conduct! — Method
Calculate conductive heat fluxes using Fourier's law
Updates array of atmos.flux_cdct at each layer of the atmosphere.
Arguments:
atmos::Atmos_tthe atmosphere struct instance to be used
AGNI.energy.convection! — Method
Calculate dry convective fluxes using mixing length theory.
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.
Uses the mixing length formulation outlined by Joyce & Tayar (2023), which was also implemented in Lee et al. (2024), and also partially outlined in the review 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
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)
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. https://doi.org/10.1029/JZ067i008p03095
The scale height is formulated as: Hp = P / (ρ g) Where ρ is obtained from the equation of state.
To account for convective stability due to compositional gradients, we can use the Ledoux criterion rather than the Schwarzschild criterion. This is described nicely in Gabriel et al. (2014), as well as Salaris & Cassisi (2017). http://dx.doi.org/10.1051/0004-6361/201423442 https://doi.org/10.1098/rsos.170192
In the ideal gas regime, the Ledoux criterion can be simply written as: ∇_ld = ∇_ad + dln(μ)/dln(P) Using Equation (10) of Gabriel+14. Taking β=Pg/P=1 means the gas pressure equals the total pressure, neglecting pressure contributions from ions and electrons.
Arguments:
atmos::Atmos_tthe atmosphere struct instance to be used.
AGNI.energy.eval_exchange_coeff — Method
Calculate turbulent kinetic energy (TKE) exchange coefficient.
Based on Monin–Obukhov similarity theory, from roughness length scale. See eq 9 in Nicholson & Benn (2009). Added small epsilon-factor to avoid function blowing-up around regime where height ≈ roughness.
Arguments:
height::Float64Height above surface [m]roughness::Float64Roughness length scale [m]
Returns:
C_d::Float64TKE exchange coefficient [dimensionless]
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_tthe atmosphere struct instance to be used.
AGNI.energy.latent! — Method
Analytical diffusion scheme for condensation and evaporation energy.
Updates fluxes. Requires chemistry.rainout_and_evaporate to be called first.
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.
If evaporation is enabled, then integrates from top downwards to determine flux from re-evaporation of droplets. Any droplets which reach the ground go towards forming an ocean.
Should ideally perform a microphysical treatment; e.g. by following this paper: https://agupubs.onlinelibrary.wiley.com/doi/epdf/10.1029/2020JE006653
Arguments:
atmos::Atmos_tthe atmosphere struct instance to be used.
AGNI.energy.make_finite! — Method
Set non-finite values in an array equal to a given fill value.
Arguments:
arrarray potentially containing non-finite valuesfillreplacement value to fill with
AGNI.energy.radtrans! — Method
Calculate radiative fluxes using the desired scheme.
Uses the configuration inside the atmos struct. Can either do LW or SW calculation, set by lw function argument.
Arguments:
atmos::Atmos_tthe atmosphere struct instance to be used.lw::Boollongwave calculation? Else: shortwavecalc_cf::Bool=falsealso calculate contribution function?
AGNI.energy.reset_fluxes! — Method
Reset energy fluxes to zero.
AGNI.energy.sensible! — Method
Calculate sensible heat flux from turbulent kinetic energy (TKE)
Updates the values of atmos.C_d and atmos.flux_sens.
Arguments:
atmos::Atmos_tthe atmosphere struct instance to be used
AGNI.energy.skin_flux — Method
Calculate flux carried by conductive skin.
Chemistry and phase-change calculations
Functions from chemistry.jl.
AGNI.chemistry — Module
This module handles chemistry, condensation, and evaporation.
Note the important distinctions between the variables which store atmospheric composition.
gas_ovmrstores VMRs inputted by the user, which are usually constant in height.gas_vmrstores the runtime gas volume mixing ratios, after all calculations are performed.
AGNI.chemistry._chem_gas! — Method
Calculate gas-phase composition at dry thermochemical equilibrium.
Uses FastChem to calculate the gas composition at each level of the atmosphere. Volatiles are converted to bulk elemental abundances, which are then provided to FastChem alongside the temperature/pressure profile. FastChem is currently called as an executable, which is not optimal.
Arguments:
atmos::Atmos_tthe atmosphere struct instance to be used.write_cfg::Boolwrite config and elements
Returns:
state::Intfastchem state (0: success, 1: criticalfail, 2: elemfail, 3: convfail, 4: bothfail)
AGNI.chemistry._sat_aloft! — Method
Handle sub/super-saturation aloft; required for latent heat flux calcaulation
Adjust gas VMRs according to saturation and cold-trap requirements.
Volatiles which are allowed to condense are rained-out at condensing levels until the gas is exactly saturated, not supersaturated. If evaporation is enabled here, it will lead to enhanced mixing ratios at deeper levels as rain is converted back into gas from a liquid state.
Any rain that reaches the surface without being evaporated leads to ocean formation. Disabling evaporation will lead to unclosed energy budget when calculation of latent heat fluxes is performed.
This function can be called after the fastchem calculation.
Arguments:
atmos::Atmos_tthe atmosphere struct instance to be used.
AGNI.chemistry._sat_surf! — Method
Handle sub/super-saturation at the surface, forming oceans.
NOTE: this function does not reset to atmos to its original state before operating.
If a condensable is supersaturated at the surface, its partial surface pressure is reduced to exactly saturation. The condensed mass is added to the ocean reservoir. If a condensable is subsaturated at the surface, its partial pressure is increased to exactly saturation, subject to the amount of ocean available. Reduces or increases the total surface pressure accordingly.
Accounts for mmw ratio when converting from partial pressures to masses: mi = (pi/g) * (mui / mutot) where mi has units [kg/m^2] and mutot is the total mixture MMW.
Arguments:
atmos::Atmos_tthe atmosphere struct instance to be used.
Returns:
any_changed::Booldid any component condense or evaporate?
AGNI.chemistry.calc_composition! — Method
Run condensation and chemistry schemes as required.
This function is designed as a wrapper for appropriately handling these three schemes together in the correct order, so that variables are appropriately updated. Steps:
- call
_sat_surf!to handle saturation at the surface, and adjust the pressure grid - call
_chem_gas!to handle gas-phase chemistry in the column - call
_sat_aloft!to handle saturation aloft, above the surface.
Arguments:
atmos::Atmos_tthe atmosphere struct instance to be used.do_surf::Booldo saturation cond/evap at surfacedo_chem::Booldo thermochemistry in the columndo_aloft::Booldo saturation cond/evap aloft
Returns:
state::Intfastchem state (0: success, 1: criticalfail, 2: elemfail, 3: convfail, 4: bothfail)
AGNI.chemistry.normalise_vmrs! — Method
Normalise gas VMRs, keeping condensates unchanged
Only acts on a single model level.
Parameters:
atmos::atmosphere.Atmos_tatmosphere structurei::Intlevel index to normalise
AGNI.chemistry.restore_composition! — Method
Reset mixing ratios, pressures, and oceans to their original values
Numerical solvers
Functions from solver.jl.
AGNI.solver.gs_search — Method
Golden section search algorithm
Minimises a function f between the bounds a and b. The function f must return a positive scalar.
Arguments:
fFunction to be minimisedaInitial bracket, lower valuebInitial bracket, upper valuiedxtolConvergence: exit when bracket is smaller than this sizeatolConvergence: absolute tolerance on minimummax_stepsMaximum number of iterationswarningsPrint warning if search does not converge beforemax_steps` are taken.
Returns:
solBest solution found by search.succDid search converge?
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_tthe atmosphere struct instance to be used.
Optional physics arguments:
sol_type::Intsolution type, 1: tmpsurf | 2: skin | 3: fluxint | 4: tgt_olrchem::Boolinclude eqm thermochemistry when solving for RCE?convect::Boolinclude convectionsens_heat::Boolinclude sensible heating at the surfaceconduct::Boolinclude conductive heat transport within the atmospherelatent::Boolinclude latent heat exchange (condensation/evaporation)advect::Boolinclude advective heat exchange (dynamics)rainout::Boolallow rainout (phase change impacts mixing ratios, not just energy fluxes)oceans::Boolcheck surface saturation (ocean formation)
Optional solver arguments:
dx_min::Float64minimum step size [K]dx_max::Float64maximum step size [K]tmp_pad::Float64padding around hard limits on temperature floor & ceiling valuesmax_steps::Intmaximum number of solver stepsmax_runtime::Float64maximum runtime in wall-clock secondsfdc::Boolfinite difference: central difference? otherwise use forward differencefdo::Intfinite difference: scheme order (2nd or 4th)method::Intnumerical method (1: Newton-Raphson, 2: Gauss-Newton, 3: Levenberg-Marquardt)easy_start::Boolimprove convergence reliability; introduce convection and latent heat graduallygrey_start::Boolimprove convergence reliability; obtain double-grey solution firstls_method::Intlinesearch algorithm (0: None, 1: golden, 2: backtracking)detect_plateau::Boolassist solver when it is stuck in a region of small dF/dTperturb_all::Boolalways recalculate entire Jacobian matrix? Otherwise updates columns only as requiredmodplot::Intiteration frequency at which to make plotssave_frames::Boolsave plotting framesmodprint::Intiteration frequency at which to print infoplot_jacobian::Boolplot jacobian too?radiative_Kzz::Boolestimate Kzz in radiative zones? Otherwise left at Kzz=0 in these regions.conv_atol::Float64convergence: absolute tolerance on per-level flux deviation [W m-2]conv_rtol::Float64convergence: relative tolerance on per-level flux deviation [dimensionless]
Returns: Nothing
AGNI.solver.solve_prescribed! — Method
Solve for global (not local) balance with a prescribed atmosphere structure.
Comparable to solve_transparent, but with an opaque prescribed atmospheric structure.
Arguments:
atmos::Atmos_tthe atmosphere struct instance to be used.sol_type::Intsolution types, same as solve_energyatm_type::Intatmosphere prescription (1: isothermal, 2: adiabat, 3: adiabat+stratosphere)conv_atol::Float64convergence: absolute tolerance on global flux [W m-2]conv_rtol::Float64convergence: relative tolerance on global flux [dimensionless]max_steps::Intmaximum number of solver stepstmp_upper::Float64upper-bound on Tsurf for golden-section search [K]
Returns:
Boolindicating success
AGNI.solver.solve_transparent! — Method
Solve for energy balance with a transparent atmosphere.
This will use an isothermal temperature profile, and only modify T_surf.
Arguments:
atmos::Atmos_tthe atmosphere struct instance to be used.sol_type::Intsolution types, same as solve_energyconv_atol::Float64convergence: absolute tolerance on global flux [W m-2]conv_rtol::Float64convergence: relative tolerance on global flux [dimensionless]max_steps::Intmaximum number of solver stepstmp_upper::Float64upper-bound on Tsurf for golden-section search [K]
Returns:
Boolindicating success
Plotting functions and utilities
Functions from plotting.jl.
AGNI.plotting.animate — Method
Convert still-frame images into an animation
Uses FFMPEG provided by Julia library.
Arguments:
output_dir::Stringfolder in which to save the animationframes_dir::Stringfolder containing the frames
Optional arguments:
output_fmt::Stringoutput file format: mp4, mov, ...frames_fmt::Stringinput frame format: png, jpg, ...duration::Float64animation duration [seconds]
AGNI.plotting.combined — Method
Combined plot used for tracking behaviour of the solver
AGNI.plotting.jacobian — Method
Plot jacobian matrix
AGNI.plotting.plot_albedo — Method
Plot spectral albedo (ratio of SWUP to SWDN)
AGNI.plotting.plot_cloud — Method
Plot the cloud mass mixing ratio and area fraction.
AGNI.plotting.plot_contfunc1 — Method
Plot contribution function at different bands.
AGNI.plotting.plot_contfunc2 — Method
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.
AGNI.plotting.plot_emission — Method
Plot emission spectrum at the TOA
AGNI.plotting.plot_fluxes — Method
Plot the fluxes at each pressure level
AGNI.plotting.plot_pt — Method
Plot the temperature-pressure profile.
AGNI.plotting.plot_radius — Method
Plot the radius vs pressure profile.
AGNI.plotting.plot_vmr — Method
Plot the VMRs of the atmosphere at each cell-centre location.
Physical properties and constants
Functions from phys.jl.
AGNI.phys._get_mmw — Method
Calculate species mean molecular weight [kg mol-1] from formula or use known value
AGNI.phys._pretty_color — Method
Generate a colour hex code from a molecular formula
AGNI.phys._pretty_name — Method
Convert formula to pretty unicode string
AGNI.phys._rho_ideal — Method
Evaluate the density of a single gas using the ideal gas EOS.
Arguments:
tmp::Float64temperature [K]prs::Float64pressure [Pa]mmw::Float64mean molecular weight [kg mol-1]
Returns:
rho::Float64mass density [kg m-3]
AGNI.phys.calc_Teq — Method
Calculate planetary equilibrium temperature.
https://en.wikipedia.org/wiki/Planetaryequilibriumtemperature?useskin=vector
Arguments:
S::Float64Bolometric instellation [W m-2]α::Float64Bond albedo
Returns:
Teq::Float64Planetary equilibrium temperature [K]
AGNI.phys.calc_Tskin — Method
Calculate planetary skin temperature.
https://en.wikipedia.org/wiki/Skintemperature(atmosphere)?useskin=vector
Arguments:
S::Float64Bolometric instellation [W m-2]α::Float64Bond albedo
Returns:
Tskin::Float64Planetary skin temperature [K]
AGNI.phys.calc_rho_gas — Method
Calculate the density of a gas using the most appropriate equation of state.
Arguments:
tmp::Float64temperature [K]prs::Float64pressure [Pa]gas::Gas_tthe gas struct to be used
Returns:
rho::Float64mass density [kg m-3]
AGNI.phys.calc_rho_mix — Method
Calculate the density of a mixture of gases using Amagat's law.
Arguments:
gas::Array{Gas_t,1}array of gasesvmr::Array{Float64,1}array of volume mixing ratiostmp::Float64temperature [K]prs::Float64pressure [Pa]
Returns:
rho::Float64mass density [kg m-3]
AGNI.phys.calc_therm_diffus — Method
Calculate thermal diffusivity.
https://en.wikipedia.org/wiki/Thermal_diffusivity?useskin=vector
Arguments:
k::Float64Thermal conductivity [W m-1 K-1]ρ::Float64Density [kg m-3]cp::Float64Specific heat capacity [J K-1 kg-1]
Returns:
α::Float64Thermal diffusivity [m2 s-1]
AGNI.phys.count_atoms — Method
Get number of atoms from formula, returning a dictionary
AGNI.phys.evaluate_planck — Method
Evaluate the Planck function at a given wavelength and temperature.
Integrated over a hemisphere.
Arguments:
wave::Float64Wavelength [nm]tmp::Float64Temperature [K]
Returns:
flx::Float64Spectral flux density [W m-2 nm-1]
AGNI.phys.get_Cp — Method
Get gas heat capacity for a given temperature.
Evaluates at 0 Celcius if gas.tmp_dep=false.
Arguments:
gas::Gas_tthe gas struct to be usedt::Float64temperature [K]
Returns:
cp::Float64heat capacity of gas [J K-1 kg-1]
AGNI.phys.get_Kc — Method
Get gas thermal conductivity at a given temperature.
This assumes that the gas is within the ideal regime.
Arguments:
gas::Gas_tthe gas struct to be usedt::Float64temperature [K]
Returns:
kc::Float64thermal conductivity [W m-1 K-1]
AGNI.phys.get_Lv — Method
Get gas enthalpy (latent heat) of phase change.
If the temperature is above the critical point, then a zero value is returned. Evaluates at 0 Celcius if gas.tmp_dep=false.
Arguments:
gas::Gas_tthe gas struct to be usedt::Float64temperature [K]
Returns:
h::Float64enthalpy of phase change [J kg-1]
AGNI.phys.get_Psat — Method
Get gas saturation pressure for a given temperature.
If the temperature is above the critical point, then a large value is returned.
Arguments:
gas::Gas_tthe gas struct to be usedt::Float64temperature [K]
Returns:
p::Float64saturation pressure [Pa]
AGNI.phys.get_Tdew — Method
Get gas dew point temperature for a given partial pressure.
If the pressure is below the critical point pressure, then T_crit is returned. This function is horrendous, and should be avoided at all costs.
Arguments:
gas::Gas_tthe gas struct to be usedp::Float64pressure [Pa]
Returns:
t::Float64dew point temperature [K]
AGNI.phys.grav_accel — Method
Calculate gravitational acceleration.
Arguments:
mass::Float64Enclosed mass [kg]radius::Float64Enclosed radius [m]
Returns:
grav::Float64Gravitational acceleration [m s-2]
AGNI.phys.is_vapour — Method
Check if pressure-temperature coordinate is within the vapour regime.
Returns true if p > psat and t < tcrit.
Arguments:
gas::Gas_tthe gas struct to be usedt::Float64temperature [K]p::Float64temperature [K]
Returns:
vapour::Boolis within vapour regime?
AGNI.phys.liquid_rho — Method
Evaluate the density of a liquid phase.
Returns BIGFLOAT density for unsupported phases, to avoid divide-by-zero error
Arguments:
name::StringName of liquid
Returns:
rho::Float64Density of liquid phase [kg m-3]
AGNI.phys.load_gas — Method
Load gas data into a new struct
AGNI.phys.same_atoms — Method
Check if two gas atom dicts are equivalent
File I/O modules
Functions from save.jl and load.jl.
AGNI.save.write_fluxes — Method
Write cell-edge energy fluxes to a CSV file
Arguments
atmos::Atmos_tAtmosphere objectfname::StringFilename to write to
AGNI.save.write_ncdf — Method
Write verbose atmosphere data to a NetCDF file
Arguments
atmos::Atmos_tAtmosphere objectfname::StringFilename to write to
AGNI.save.write_profile — Method
Write {Pressure, Temperature, Radius} profile to a CSV file
Arguments
atmos::Atmos_tAtmosphere objectfname::StringFilename to write to
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_tthe atmosphere struct instance to be used.fname::Stringpath to the NetCDF file
Returns:
ok::Boolloaded ok?