pycontrails.models.cocip.contrail_properties

Contrail Property Calculations.

Functions

contrail_edges(lon, lat, sin_a, cos_a, width)

Calculate the longitude and latitude of the contrail edges to account for contrail spreading.

contrail_optical_depth(r_ice_vol, ...)

Calculate the contrail optical depth for each waypoint.

contrail_persistent(latitude, altitude, ...)

Determine surviving contrail segments after time integration step.

contrail_vertices(lon, lat, sin_a, cos_a, ...)

Calculate the longitude and latitude of the contrail vertices.

energy_forcing(rf_net_t1, rf_net_t2, ...)

Calculate the contrail energy forcing over time step.

horizontal_diffusivity(ds_dz, depth)

Calculate contrail horizontal diffusivity.

ice_particle_activation_rate(...)

Calculate the activation rate of black carbon particles to contrail ice crystals.

ice_particle_mass(r_ice_vol)

Calculate the contrail ice particle mass.

ice_particle_number(nvpm_ei_n, fuel_dist, ...)

Calculate the initial number of ice particles per distance after the wake vortex phase.

ice_particle_number_per_mass_of_air(...)

Calculate the number of contrail ice particles per mass of air.

ice_particle_number_per_volume_of_plume(...)

Calculate the number of contrail ice particles per volume of plume (n_ice_per_vol).

ice_particle_survival_fraction(iwc, iwc_1)

Estimate the fraction of contrail ice particle number that survive the wake vortex phase.

ice_particle_terminal_fall_speed(...)

Calculate the terminal fall speed of contrail ice particles.

ice_particle_volume_mean_radius(iwc, ...)

Calculate the ice particle volume mean radius.

initial_iwc(air_temperature, ...)

Estimate the initial contrail ice water content (iwc) before the wake vortex phase.

initial_persistent(iwc_1, rhi_1)

Determine if waypoints have persistent contrails.

iwc_adiabatic_heating(air_temperature, ...)

Calculate the change in ice water content due to adiabatic heating from the wake vortex phase.

iwc_post_wake_vortex(iwc, iwc_ad)

Calculate the ice water content after the wake vortex phase (iwc_1).

light_wave_phase_delay(r_ice_vol)

Calculate the phase delay of the light wave passing through the contrail ice particle.

mean_energy_flux_per_m(rad_flux_per_m, dt)

Calculate the mean energy flux per length of contrail on segment following waypoint.

mean_radiative_flux_per_m(rf_net_t1, ...)

Calculate the mean radiative flux per length of contrail between two time steps.

new_contrail_dimensions(sigma_yy_t2, sigma_zz_t2)

Calculate the new contrail width and depth.

new_effective_area_from_sigma(sigma_yy, ...)

Calculate effective cross-sectional area of contrail plume (area_eff) from sigma parameters.

new_ice_particle_number(n_ice_per_m_t1, ...)

Calculate the number of ice particles per distance at the end of the time step.

new_ice_water_content(iwc_t1, q_t1, q_t2, ...)

Calculate the new contrail ice water content after the time integration step (iwc_t2).

particle_losses_aggregation(r_ice_vol, ...)

Calculate the rate of contrail ice particle losses due to sedimentation-induced aggregation.

particle_losses_turbulence(width, depth, ...)

Calculate the rate of contrail ice particle losses due to plume-internal turbulence.

plume_effective_cross_sectional_area(width, ...)

Calculate the effective cross-sectional area of the contrail plume (area_eff).

plume_effective_depth(width, area_eff)

Calculate the effective depth of the contrail plume (depth_eff).

plume_mass_per_distance(area_eff, rho_air)

Calculate the contrail plume mass per unit length.

plume_temporal_evolution(width_t1, depth_t1, ...)

Calculate the temporal evolution of the contrail plume parameters.

q_exhaust(air_temperature, air_pressure, ...)

Calculate the specific humidity released by water vapor from aircraft emissions.

scattering_extinction_efficiency(r_ice_vol)

Calculate the scattering extinction efficiency (q_ext) based on Mie-theory.

segment_length_ratio(seg_length_t1, ...)

Calculate the ratio of contrail segment length pre-advection to post-advection.

temperature_adiabatic_heating(...)

Calculate the ambient air temperature for each waypoint after the wake vortex phase.

vertical_diffusivity(air_pressure, ...)

Calculate contrail vertical diffusivity.

pycontrails.models.cocip.contrail_properties.contrail_edges(lon, lat, sin_a, cos_a, width)

Calculate the longitude and latitude of the contrail edges to account for contrail spreading.

(lon_edge_l, lat_edge_l) x———————

(Contrail midpoint: lon, lat) X===================== ->

(lon_edge_r, lat_edge_r) x———————

Parameters:
  • lon (npt.NDArray[np.float64]) – longitude of contrail waypoint, degrees

  • lat (npt.NDArray[np.float64]) – latitude of contrail waypoint, degrees

  • sin_a (npt.NDArray[np.float64]) – sin(a), where a is the angle between the plume and the longitudinal axis

  • cos_a (npt.NDArray[np.float64]) – cos(a), where a is the angle between the plume and the longitudinal axis

  • width (npt.NDArray[np.float64]) – contrail width at each waypoint, [\(m\)]

Returns:

tuple[npt.NDArray[np.float64], npt.NDArray[np.float64], npt.NDArray[np.float64], npt.NDArray[np.float64]] – (lon_edge_l, lat_edge_l, lon_edge_r, lat_edge_r), longitudes and latitudes at the left and right edges of the contrail, degrees

pycontrails.models.cocip.contrail_properties.contrail_optical_depth(r_ice_vol, n_ice_per_m, width)

Calculate the contrail optical depth for each waypoint.

Parameters:
  • r_ice_vol (npt.NDArray[np.float64]) – ice particle volume mean radius, [\(m\)]

  • n_ice_per_m (npt.NDArray[np.float64]) – Number of contrail ice particles per distance, [\(m^{-1}\)]

  • width (npt.NDArray[np.float64]) – Contrail width, [\(m\)]

Returns:

npt.NDArray[np.float64] – Contrail optical depth

pycontrails.models.cocip.contrail_properties.contrail_persistent(latitude, altitude, segment_length, age, tau_contrail, n_ice_per_m3, params)

Determine surviving contrail segments after time integration step.

A contrail waypoint reaches its end of life if any of the following conditions hold:

  1. Contrail age exceeds max_age

  2. Contrail optical depth lies outside of interval [min_tau, max_tau]

  3. Ice number density lies outside of interval [min_n_ice_per_m3, max_n_ice_per_m3]

  4. Altitude lies outside of the interval [min_altitude_m, max_altitude_m]

  5. Segment length exceeds max_seg_length_m

  6. Latitude values are within 1 degree of the north or south pole

This function warns if all values in the tau_contrail array are nan.

Changed in version 0.25.10: Extreme values of latitude (ie, close to the north or south pole) now create an end of life condition. This check helps address issues related to divergence in the polar regions. (With large enough integration time delta, it is still possible for post-advection latitude values to lie outside of [-90, 90], but this is no longer possible with typical parameters and wind speeds.)

Parameters:
  • latitude (npt.NDArray[np.float64]) – Contrail latitude, [\(\deg\)]

  • altitude (npt.NDArray[np.float64]) – Contrail altitude, [\(m\)]

  • segment_length (npt.NDArray[np.float64]) – Contrail segment length, [\(m\)]

  • age (npt.NDArray[np.timedelta64]) – Contrail age

  • tau_contrail (npt.NDArray[np.float64]) – Contrail optical depth

  • n_ice_per_m3 (npt.NDArray[np.float64]) – Contrail ice particle number per volume of air, [\(# m^{-3}\)]

  • params (dict[str, Any]) – Dictionary of CocipParams parameters determining the conditions for end of contrail life.

Returns:

npt.NDArray[np.bool_] – Boolean array indicating surviving contrails. Persisting contrails will be marked as True.

pycontrails.models.cocip.contrail_properties.contrail_vertices(lon, lat, sin_a, cos_a, width, segment_length)

Calculate the longitude and latitude of the contrail vertices.

This is equivalent to running contrail_edges() at each contrail waypoint and associating the next continuous waypoint with the previous. This method is helpful when you want to treat each contrail waypoint independently.

(lon_1, lat_1) x——————–x (lon_4, lat_4)

(Contrail waypoint: lon, lat) X==================== ->

(lon_2, lat_2) x——————–x (lon_3, lat_3)

Parameters:
  • lon (npt.NDArray[np.float64]) – longitude of contrail waypoint, degrees

  • lat (npt.NDArray[np.float64]) – latitude of contrail waypoint, degrees

  • sin_a (npt.NDArray[np.float64]) – sin(a), where a is the angle between the plume and the longitudinal axis

  • cos_a (npt.NDArray[np.float64]) – cos(a), where a is the angle between the plume and the longitudinal axis

  • width (npt.NDArray[np.float64]) – contrail width at each waypoint, [\(m\)]

  • segment_length (npt.NDArray[np.float64]) – contrail length at each waypoint, [\(m\)]

Returns:

tuple[npt.NDArray[np.float64], npt.NDArray[np.float64], npt.NDArray[np.float64], npt.NDArray[np.float64]] – (lon_1, lat_1, lon_2, lat_2, lon_3, lat_3, lon_4, lat_4) degrees

pycontrails.models.cocip.contrail_properties.energy_forcing(rf_net_t1, rf_net_t2, width_t1, width_t2, seg_length_t2, dt)

Calculate the contrail energy forcing over time step.

The contrail energy forcing is calculated as the local contrail net radiative forcing (RF’, change in energy flux per contrail area) multiplied by its width and integrated over its length and lifetime.

Parameters:
  • rf_net_t1 (npt.NDArray[np.float64]) – local contrail net radiative forcing at the start of the time step, [\(W m^{-2}\)]

  • rf_net_t2 (npt.NDArray[np.float64]) – local contrail net radiative forcing at the end of the time step, [\(W m^{-2}\)]

  • width_t1 (npt.NDArray[np.float64]) – contrail width at the start of the time step, [\(m\)]

  • width_t2 (npt.NDArray[np.float64]) – contrail width at the end of the time step, [\(m\)]

  • seg_length_t2 (npt.NDArray[np.float64] | float) – Segment length of contrail waypoint at the end of the time step, [\(m\)]

  • dt (npt.NDArray[np.timedelta64] | np.timedelta64) – integrate contrails with time steps of dt, [\(s\)]

Returns:

npt.NDArray[np.float64] – Contrail energy forcing over time step dt, [\(J\)].

pycontrails.models.cocip.contrail_properties.horizontal_diffusivity(ds_dz, depth)

Calculate contrail horizontal diffusivity.

Parameters:
  • ds_dz (npt.NDArray[np.float64]) – Total wind shear (eastward and northward winds) with respect to altitude (dz), [\(m s^{-1} / Pa\)]

  • depth (npt.NDArray[np.float64]) – Contrail depth at each waypoint, [\(m\)]

Returns:

npt.NDArray[np.float64] – horizontal diffusivity, [\(m^{2} s^{-1}\)]

References

Notes

Accounts for the turbulence-induced diffusive contrail spreading in the horizontal direction.

pycontrails.models.cocip.contrail_properties.ice_particle_activation_rate(air_temperature, T_crit_sac)

Calculate the activation rate of black carbon particles to contrail ice crystals.

The activation rate is calculated as a function of the difference between the ambient temperature and the Schmidt-Appleman threshold temperature T_crit_sac.

Parameters:
  • air_temperature (npt.NDArray[np.float64]) – ambient temperature at each waypoint before wake_wortex, [\(K\)]

  • T_crit_sac (npt.NDArray[np.float64]) – estimated Schmidt-Appleman temperature threshold for contrail formation, [\(K\)]

Returns:

npt.NDArray[np.float64] – Proportion of black carbon particles that activates to contrail ice parties.

Notes

The equation is not published but based on the raw data from [Bräuer et al., 2021].

References

pycontrails.models.cocip.contrail_properties.ice_particle_mass(r_ice_vol)

Calculate the contrail ice particle mass.

It is calculated by multiplying the mean ice particle volume with the density of ice

Parameters:

r_ice_vol (npt.NDArray[np.float64]) – Ice particle volume mean radius, [\(m\)]

Returns:

npt.NDArray[np.float64] – Mean contrail ice particle mass, [\(kg\)]

pycontrails.models.cocip.contrail_properties.ice_particle_number(nvpm_ei_n, fuel_dist, f_surv, air_temperature, T_crit_sac, min_ice_particle_number_nvpm_ei_n)

Calculate the initial number of ice particles per distance after the wake vortex phase.

The initial number of ice particle per distance is calculated from the black carbon number emissions index nvpm_ei_n and fuel burn per distance fuel_dist. Note that a lower bound for nvpm_ei_n is set at 1e13 \(kg^{-1}\) to account for the activation of ambient aerosol particles and organic volatile particles.

Parameters:
  • nvpm_ei_n (npt.NDArray[np.float64]) – black carbon number emissions index, [\(kg^{-1}\)]

  • fuel_dist (npt.NDArray[np.float64]) – fuel consumption of the flight segment per distance traveled, [\(kg m^{-1}\)]

  • f_surv (npt.NDArray[np.float64]) – Fraction of contrail ice particle number that survive the wake vortex phase.

  • air_temperature (npt.NDArray[np.float64]) – ambient temperature for each waypoint, [\(K\)]

  • T_crit_sac (npt.NDArray[np.float64]) – estimated Schmidt-Appleman temperature threshold for contrail formation, [\(K\)]

  • min_ice_particle_number_nvpm_ei_n (float) – lower bound for nvpm_ei_n to account for ambient aerosol particles for newer engines [\(kg^{-1}\)]

Returns:

npt.NDArray[np.float64] – initial number of ice particles per distance after the wake vortex phase, [\(# m^{-1}\)]

pycontrails.models.cocip.contrail_properties.ice_particle_number_per_mass_of_air(n_ice_per_vol, rho_air)

Calculate the number of contrail ice particles per mass of air.

Parameters:
  • n_ice_per_vol (npt.NDArray[np.float64]) – number of ice particles per volume of contrail plume at time t, [\(# m^{-3}\)]

  • rho_air (npt.NDArray[np.float64]) – density of air for each waypoint, [\(kg m^{-3}\)]

Returns:

npt.NDArray[np.float64] – number of ice particles per mass of air at time t, [\(# kg^{-1}\)]

pycontrails.models.cocip.contrail_properties.ice_particle_number_per_volume_of_plume(n_ice_per_m, area_eff)

Calculate the number of contrail ice particles per volume of plume (n_ice_per_vol).

Parameters:
  • n_ice_per_m (npt.NDArray[np.float64]) – number of ice particles per distance at time t, [\(m^{-1}\)]

  • area_eff (npt.NDArray[np.float64]) – effective cross-sectional area of the contrail plume, [\(m^{2}\)]

Returns:

npt.NDArray[np.float64] – number of ice particles per volume of contrail plume at time t, [\(# m^{-3}\)]

pycontrails.models.cocip.contrail_properties.ice_particle_survival_fraction(iwc, iwc_1)

Estimate the fraction of contrail ice particle number that survive the wake vortex phase.

CoCiP assumes that this fraction is proportional to the change in ice water content (iwc_1 - iwc) before and after the wake vortex phase.

Parameters:
  • iwc (npt.NDArray[np.float64]) – initial ice water content at each waypoint before the wake vortex phase, [\(kg_{H_{2}O}/kg_{air}\)]

  • iwc_1 (npt.NDArray[np.float64]) – ice water content after the wake vortex phase, [\(kg_{H_{2}O}/kg_{air}\)]

Returns:

npt.NDArray[np.float64] – Fraction of contrail ice particle number that survive the wake vortex phase.

pycontrails.models.cocip.contrail_properties.ice_particle_terminal_fall_speed(air_pressure, air_temperature, r_ice_vol)

Calculate the terminal fall speed of contrail ice particles.

v_t is calculated based on a parametric model from [Spichtinger and Gierens, 2009], using inputs of pressure level, ambient temperature and the ice particle volume mean radius. See Table 2 for the model parameters.

Parameters:
  • air_pressure (npt.NDArray[np.float64]) – Pressure altitude at each waypoint, [\(Pa\)]

  • air_temperature (npt.NDArray[np.float64]) – Ambient temperature for each waypoint, [\(K\)]

  • r_ice_vol (npt.NDArray[np.float64]) – Ice particle volume mean radius, [\(m\)]

Returns:

npt.NDArray[np.float64] – Terminal fall speed of contrail ice particles, [\(m s^{-1}\)]

References

pycontrails.models.cocip.contrail_properties.ice_particle_volume_mean_radius(iwc, n_ice_per_kg_air)

Calculate the ice particle volume mean radius.

Parameters:
  • iwc (npt.NDArray[np.float64]) – contrail ice water content, i.e., contrail ice mass per kg of air, [\(kg_{H_{2}O}/kg_{air}\)]

  • n_ice_per_kg_air (npt.NDArray[np.float64]) – number of ice particles per mass of air, [\(# kg^{-1}\)]

Returns:

npt.NDArray[np.float64] – ice particle volume mean radius, [\(m\)]

Notes

r_ice_vol is the mean radius of a sphere that has the same volume as the contrail ice particle.

r_ice_vol calculated by dividing the total volume of contrail ice particle per kg of air (total_ice_volume, \(m**3/kg-air\)) with the number of contrail ice particles per kg of air (n_ice_per_kg_air, \(#/kg-air\)).

pycontrails.models.cocip.contrail_properties.initial_iwc(air_temperature, specific_humidity, air_pressure, fuel_dist, width, depth, ei_h2o)

Estimate the initial contrail ice water content (iwc) before the wake vortex phase.

Note that the ice water content is replaced by zero if it is negative (dry air), and this will end the contrail life-cycle in subsequent steps.

Parameters:
  • air_temperature (npt.NDArray[np.float64]) – ambient temperature for each waypoint, [\(K\)]

  • specific_humidity (npt.NDArray[np.float64]) – ambient specific humidity for each waypoint, [\(kg_{H_{2}O}/kg_{air}\)]

  • air_pressure (npt.NDArray[np.float64]) – initial pressure altitude at each waypoint, before the wake vortex phase, [\(Pa\)]

  • fuel_dist (npt.NDArray[np.float64]) – fuel consumption of the flight segment per distance traveled, [\(kg m^{-1}\)]

  • width (npt.NDArray[np.float64]) – initial contrail width, [\(m\)]

  • depth (npt.NDArray[np.float64]) – initial contrail depth, [\(m\)]

  • ei_h2o (float) – water vapor emissions index of fuel, [\(kg_{H_{2}O} \ kg_{fuel}^{-1}\)]

Returns:

npt.NDArray[np.float64] – Initial contrail ice water content (iwc) at the original waypoint before the wake vortex phase, [\(kg_{H_{2}O}/kg_{air}\)]. Returns zero if iwc is is negative (dry air).

pycontrails.models.cocip.contrail_properties.initial_persistent(iwc_1, rhi_1)

Determine if waypoints have persistent contrails.

Conditions for persistent initial_contrails:

  1. ice water content at level 1: 1e-12 < iwc < 1e10

  2. rhi at level 1: 0 < rhi < 1e10

Changed in version 0.25.1: Returned array now has floating dtype. This is consistent with other filtering steps in the CoCiP model (ie, sac).

Parameters:
  • iwc_1 (npt.NDArray[np.float64]) – ice water content after the wake vortex phase, [\(kg_{H_{2}O}/kg_{air}\)]

  • rhi_1 (npt.NDArray[np.float64]) – relative humidity with respect to ice after the wake vortex phase

Returns:

npt.NDArray[np.float64] – Mask of waypoints with persistent contrails. Waypoints with persistent contrails will have value 1.

Notes

The RHi at level 1 does not need to be above 100% if the iwc > 0 kg/kg. If iwc > 0 and RHi < 100%, the contrail lifetime will not end immediately as the ice particles will gradually evaporate with the rate depending on the background RHi.

pycontrails.models.cocip.contrail_properties.iwc_adiabatic_heating(air_temperature, air_pressure, air_pressure_1)

Calculate the change in ice water content due to adiabatic heating from the wake vortex phase.

Parameters:
  • air_temperature (npt.NDArray[np.float64]) – ambient temperature for each waypoint, [\(K\)]

  • air_pressure (npt.NDArray[np.float64]) – initial pressure altitude at each waypoint, before the wake vortex phase, [\(Pa\)]

  • air_pressure_1 (npt.NDArray[np.float64]) – pressure altitude at each waypoint, after the wake vortex phase, [\(Pa\)]

Returns:

npt.NDArray[np.float64] – Change in ice water content due to adiabatic heating from the wake vortex phase, [\(kg_{H_{2}O}/kg_{air}\)]

pycontrails.models.cocip.contrail_properties.iwc_post_wake_vortex(iwc, iwc_ad)

Calculate the ice water content after the wake vortex phase (iwc_1).

iwc_1 is calculated by subtracting the initial iwc before the wake vortex phase (iwc) by the change in iwc from adiabatic heating experienced during the wake vortex phase.

Note that the iwc is replaced by zero if it is negative (dry air), and this will end the contrail lifecycle in subsequent steps.

Parameters:
  • iwc (npt.NDArray[np.float64]) – initial ice water content at each waypoint before the wake vortex phase, [\(kg_{H_{2}O}/kg_{air}\)]

  • iwc_ad (npt.NDArray[np.float64]) – change in iwc from adiabatic heating during the wake vortex phase, [\(kg_{H_{2}O}/kg_{air}\)]

Returns:

npt.NDArray[np.float64] – ice water content after the wake vortex phase, iwc_1, [\(kg_{H_{2}O}/kg_{air}\)]

Notes

Level 1, see Figure 1 of [Schumann, 2012]

References

pycontrails.models.cocip.contrail_properties.light_wave_phase_delay(r_ice_vol)

Calculate the phase delay of the light wave passing through the contrail ice particle.

Parameters:

r_ice_vol (npt.NDArray[np.float64]) – ice particle volume mean radius, [\(m\)]

Returns:

npt.NDArray[np.float64] – phase delay of the light wave passing through the contrail ice particle

References

pycontrails.models.cocip.contrail_properties.mean_energy_flux_per_m(rad_flux_per_m, dt)

Calculate the mean energy flux per length of contrail on segment following waypoint.

Parameters:
  • rad_flux_per_m (npt.NDArray[np.float64]) – Mean radiative flux between time steps for waypoint, [\(W m^{-1}\)]. See mean_radiative_flux_per_m().

  • dt (npt.NDArray[np.timedelta64]) – timedelta of integration timestep for each waypoint.

Returns:

npt.NDArray[np.float64] – Mean energy flux per length of contrail after waypoint, [\(J m^{-1}\)]

Notes

Implementation differs from original fortran in two ways:

  • Discontinuity is no longer set to 0 (this occurs directly in model Cocip)

  • Instead of taking an average of the previous and following segments, energy flux is only calculated for the following segment.

pycontrails.models.cocip.contrail_properties.mean_radiative_flux_per_m(rf_net_t1, rf_net_t2, width_t1, width_t2)

Calculate the mean radiative flux per length of contrail between two time steps.

Parameters:
  • rf_net_t1 (npt.NDArray[np.float64]) – local contrail net radiative forcing at the start of the time step, [\(W m^{-2}\)]

  • rf_net_t2 (npt.NDArray[np.float64]) – local contrail net radiative forcing at the end of the time step, [\(W m^{-2}\)]

  • width_t1 (npt.NDArray[np.float64]) – contrail width at the start of the time step, [\(m\)]

  • width_t2 (npt.NDArray[np.float64]) – contrail width at the end of the time step, [\(m\)]

Returns:

npt.NDArray[np.float64] – Mean radiative flux between time steps, [\(W m^{-1}\)]

pycontrails.models.cocip.contrail_properties.new_contrail_dimensions(sigma_yy_t2, sigma_zz_t2)

Calculate the new contrail width and depth.

Parameters:
  • sigma_yy_t2 (npt.NDArray[np.float64]) – element yy, covariance matrix of the Gaussian concentration field, Eq. (6) of Schumann (2012)

  • sigma_zz_t2 (npt.NDArray[np.float64]) – element zz, covariance matrix of the Gaussian concentration field, Eq. (6) of Schumann (2012)

Returns:

  • width_t2 (npt.NDArray[np.float64]) – Contrail width at the end of the time step, [\(m\)]

  • depth_t2 (npt.NDArray[np.float64]) – Contrail depth at the end of the time step, [\(m\)]

pycontrails.models.cocip.contrail_properties.new_effective_area_from_sigma(sigma_yy, sigma_zz, sigma_yz)

Calculate effective cross-sectional area of contrail plume (area_eff) from sigma parameters.

This method calculates the same output as :func`plume_effective_cross_sectional_area`, but calculated with different input parameters.

Parameters:
  • sigma_yy (npt.NDArray[np.float64]) – element yy, covariance matrix of the Gaussian concentration field, Eq. (6) of Schumann (2012)

  • sigma_zz (npt.NDArray[np.float64]) – element zz, covariance matrix of the Gaussian concentration field, Eq. (6) of Schumann (2012)

  • sigma_yz (npt.NDArray[np.float64] | float) – element yz, covariance matrix of the Gaussian concentration field, Eq. (6) of Schumann (2012)

Returns:

npt.NDArray[np.float64] – Effective cross-sectional area of the contrail plume (area_eff)

pycontrails.models.cocip.contrail_properties.new_ice_particle_number(n_ice_per_m_t1, dn_dt_agg, dn_dt_turb, seg_ratio, dt)

Calculate the number of ice particles per distance at the end of the time step.

Parameters:
  • n_ice_per_m_t1 (npt.NDArray[np.float64]) – number of contrail ice particles per distance at the start of the time step, [\(m^{-1}\)]

  • dn_dt_agg (npt.NDArray[np.float64]) – rate of ice particle losses due to sedimentation-induced aggregation, [\(# s^{-1}\)]

  • dn_dt_turb (npt.NDArray[np.float64]) – rate of contrail ice particle losses due to plume-internal turbulence, [\(# s^{-1}\)]

  • seg_ratio (npt.NDArray[np.float64] | float) – Segment length ratio before and after it is advected to the new location.

  • dt (npt.NDArray[np.timedelta64] | np.timedelta64) – integrate contrails with time steps of dt, [\(s\)]

Returns:

npt.NDArray[np.float64] – number of ice particles per distance at the end of the time step, [\(m^{-1}\)]

pycontrails.models.cocip.contrail_properties.new_ice_water_content(iwc_t1, q_t1, q_t2, q_sat_t1, q_sat_t2, mass_plume_t1, mass_plume_t2)

Calculate the new contrail ice water content after the time integration step (iwc_t2).

Parameters:
  • iwc_t1 (npt.NDArray[np.float64]) – contrail ice water content, i.e., contrail ice mass per kg of air, at the start of the time step, [\(kg_{H_{2}O}/kg_{air}\)]

  • q_t1 (npt.NDArray[np.float64]) – specific humidity for each waypoint at the start of the time step, [\(kg_{H_{2}O}/kg_{air}\)]

  • q_t2 (npt.NDArray[np.float64]) – specific humidity for each waypoint at the end of the time step, [\(kg_{H_{2}O}/kg_{air}\)]

  • q_sat_t1 (npt.NDArray[np.float64]) – saturation humidity for each waypoint at the start of the time step, [\(kg_{H_{2}O}/kg_{air}\)]

  • q_sat_t2 (npt.NDArray[np.float64]) – saturation humidity for each waypoint at the end of the time step, [\(kg_{H_{2}O}/kg_{air}\)]

  • mass_plume_t1 (npt.NDArray[np.float64]) – contrail plume mass per unit length at the start of the time step, [\(kg_{air} m^{-1}\)]

  • mass_plume_t2 (npt.NDArray[np.float64]) – contrail plume mass per unit length at the end of the time step, [\(kg_{air} m^{-1}\)]

Returns:

npt.NDArray[np.float64] – Contrail ice water content at the end of the time step, [\(kg_{ice} kg_{air}^{-1}\)]

Notes

  1. The ice water content is fully conservative.

  2. mass_h2o_t2: the total H2O mass (ice + vapour) per unit of contrail plume [Units of kg-H2O/m]

  3. q_sat is used to calculate mass_h2o because air inside the contrail is assumed to be ice saturated.

  4. (mass_plume_t2 - mass_plume) * q_mean: contrail absorbs (releases) H2O from (to) surrounding air.

  5. iwc_t2 = mass_h2o_t2 / mass_plume_t2 - q_sat_t2: H2O in the gas phase is removed (- q_sat_t2).

pycontrails.models.cocip.contrail_properties.particle_losses_aggregation(r_ice_vol, terminal_fall_speed, area_eff, agg_efficiency=1.0)

Calculate the rate of contrail ice particle losses due to sedimentation-induced aggregation.

Parameters:
  • r_ice_vol (npt.NDArray[np.float64]) – Ice particle volume mean radius, [\(m\)]

  • terminal_fall_speed (npt.NDArray[np.float64]) – Terminal fall speed of contrail ice particles, [\(m s^{-1}\)]

  • area_eff (npt.NDArray[np.float64]) – Effective cross-sectional area of the contrail plume, [\(m^{2}\)]

  • agg_efficiency (float, optional) – Aggregation efficiency

Returns:

npt.NDArray[np.float64] – Rate of contrail ice particle losses due to sedimentation-induced aggregation, [\(# s^{-1}\)]

Notes

The aggregation efficiency (agg_efficiency = 1) was calibrated based on the observed lifetime and optical properties from the Contrail Library (COLI) database ([Schumann et al., 2017]).

References

pycontrails.models.cocip.contrail_properties.particle_losses_turbulence(width, depth, depth_eff, diffuse_h, diffuse_v, turb_efficiency=0.1)

Calculate the rate of contrail ice particle losses due to plume-internal turbulence.

Parameters:
  • width (npt.NDArray[np.float64]) – Contrail width at each waypoint, [\(m\)]

  • depth (npt.NDArray[np.float64]) – Contrail depth at each waypoint, [\(m\)]

  • depth_eff (npt.NDArray[np.float64]) – Effective depth of the contrail plume, [\(m\)]

  • diffuse_h (npt.NDArray[np.float64]) – Horizontal diffusivity, [\(m^{2} s^{-1}\)]

  • diffuse_v (npt.NDArray[np.float64]) – Vertical diffusivity, [\(m^{2} s^{-1}\)]

  • turb_efficiency (float, optional) – Turbulence sublimation efficiency

Returns:

npt.NDArray[np.float64] – Rate of contrail ice particle losses due to plume-internal turbulence, [\(# s^{-1}\)]

Notes

The turbulence sublimation efficiency (turb_efficiency = 0.1) was calibrated based on the observed lifetime and optical properties from the Contrail Library (COLI) database ([Schumann et al., 2017]).

References

pycontrails.models.cocip.contrail_properties.plume_effective_cross_sectional_area(width, depth, sigma_yz)

Calculate the effective cross-sectional area of the contrail plume (area_eff).

sigma_yy, sigma_zz and sigma_yz are the parameters governing the contrail plume’s temporal evolution.

Parameters:
  • width (npt.NDArray[np.float64]) – contrail width at each waypoint, [\(m\)]

  • depth (npt.NDArray[np.float64]) – contrail depth at each waypoint, [\(m\)]

  • sigma_yz (npt.NDArray[np.float64] | float) – temporal evolution of the contrail plume parameters

Returns:

npt.NDArray[np.float64] – effective cross-sectional area of the contrail plume, [\(m^{2}\)]

pycontrails.models.cocip.contrail_properties.plume_effective_depth(width, area_eff)

Calculate the effective depth of the contrail plume (depth_eff).

depth_eff is calculated from the effective cross-sectional area (area_eff) and the contrail width.

Parameters:
  • width (npt.NDArray[np.float64]) – contrail width at each waypoint, [\(m\)]

  • area_eff (npt.NDArray[np.float64]) – effective cross-sectional area of the contrail plume, [\(m^{2}\)]

Returns:

npt.NDArray[np.float64] – effective depth of the contrail plume, [\(m\)]

pycontrails.models.cocip.contrail_properties.plume_mass_per_distance(area_eff, rho_air)

Calculate the contrail plume mass per unit length.

Parameters:
  • area_eff (npt.NDArray[np.float64]) – effective cross-sectional area of the contrail plume, [\(m^{2}\)]

  • rho_air (npt.NDArray[np.float64]) – density of air for each waypoint, [\(kg m^{-3}\)]

Returns:

npt.NDArray[np.float64] – contrail plume mass per unit length, [\(kg m^{-1}\)]

pycontrails.models.cocip.contrail_properties.plume_temporal_evolution(width_t1, depth_t1, sigma_yz_t1, dsn_dz_t1, diffuse_h_t1, diffuse_v_t1, seg_ratio, dt, max_depth)

Calculate the temporal evolution of the contrail plume parameters.

Refer to equation (6) of Schumann (2012). See also equations (29), (30), and (31).

Parameters:
  • width_t1 (npt.NDArray[np.float64]) – contrail width at the start of the time step, [\(m\)]

  • depth_t1 (npt.NDArray[np.float64]) – contrail depth at the start of the time step, [\(m\)]

  • sigma_yz_t1 (npt.NDArray[np.float64]) – sigma_yz governs the contrail plume’s temporal evolution at the start of the time step

  • dsn_dz_t1 (npt.NDArray[np.float64]) – vertical gradient of the horizontal velocity (wind shear) normal to the contrail axis at the start of the time step, [\(m s^{-1} / Pa\)]:

    X-----------------------X               X
                 ^                           |
                 | (dsn_dz)                  |  <-- (dsn_dz)
                 |                           |
                                             X
    
  • diffuse_h_t1 (npt.NDArray[np.float64]) – horizontal diffusivity at the start of the time step, [\(m^{2} s^{-1}\)]

  • diffuse_v_t1 (npt.NDArray[np.float64]) – vertical diffusivity at the start of the time step, [\(m^{2} s^{-1}\)]

  • seg_ratio (npt.NDArray[np.float64] | float) – Segment length ratio before and after it is advected to the new location. See segment_length_ratio().

  • dt (npt.NDArray[np.timedelta64] | np.timedelta64) – integrate contrails with time steps of dt, [\(s\)]

  • max_depth (float | None) – Constrain maximum plume depth to prevent unrealistic values, [\(m\)]. If None is passed, the maximum plume depth is not constrained.

Returns:

  • sigma_yy_t2 (npt.NDArray[np.float64]) – The yy component of convariance matrix, [\(m^{2}\)]

  • sigma_zz_t2 (npt.NDArray[np.float64]) – The zz component of convariance matrix, [\(m^{2}\)]

  • sigma_yz_t2 (npt.NDArray[np.float64]) – The yz component of convariance matrix, [\(m^{2}\)]

pycontrails.models.cocip.contrail_properties.q_exhaust(air_temperature, air_pressure, fuel_dist, width, depth, ei_h2o)

Calculate the specific humidity released by water vapor from aircraft emissions.

Parameters:
  • air_temperature (npt.NDArray[np.float64]) – ambient temperature for each waypoint, [\(K\)]

  • air_pressure (npt.NDArray[np.float64]) – initial pressure altitude at each waypoint, before the wake vortex phase, [\(Pa\)]

  • fuel_dist (npt.NDArray[np.float64]) – fuel consumption of the flight segment per distance travelled, [\(kg m^{-1}\)]

  • width (npt.NDArray[np.float64]) – initial contrail width, [\(m\)]

  • depth (npt.NDArray[np.float64]) – initial contrail depth, [\(m\)]

  • ei_h2o (float) – water vapor emissions index of fuel, [\(kg_{H_{2}O} \ kg_{fuel}^{-1}\)]

Returns:

npt.NDArray[np.float64] – Humidity released by water vapour from aircraft emissions, [\(kg_{H_{2}O}/kg_{air}\)]

pycontrails.models.cocip.contrail_properties.scattering_extinction_efficiency(r_ice_vol)

Calculate the scattering extinction efficiency (q_ext) based on Mie-theory.

Parameters:

r_ice_vol (npt.NDArray[np.float64]) – ice particle volume mean radius, [\(m\)]

Returns:

npt.NDArray[np.float64] – scattering extinction efficiency

References

pycontrails.models.cocip.contrail_properties.segment_length_ratio(seg_length_t1, seg_length_t2)

Calculate the ratio of contrail segment length pre-advection to post-advection.

Parameters:
  • seg_length_t1 (npt.NDArray[np.float64]) – Segment length of contrail waypoint at the start of the time step, [\(m\)]

  • seg_length_t2 (npt.NDArray[np.float64]) – Segment length of contrail waypoint after time step and advection, [\(m\)]

Returns:

npt.NDArray[np.float64] – Ratio of segment length before advection to segment length after advection.

Notes

This implementation differs from the original fortran implementation. Instead of taking a geometric mean between the previous and following segments, a simple ratio is computed.

For terminal waypoints along a flight trajectory, the associated segment length is 0. In this case, the segment ratio is set to 1 (the naive ratio 0 / 0 is undefined). According to CoCiP conventions, terminus waypoints are “discontinuous” within the flight trajectory, and will not contribute to contrail calculations.

More broadly, any undefined (nan values, or division by 0) segment ratio is set to 1. This convention ensures that the contrail calculations are not affected by undefined segment-based properties.

Presently, the output of this function is only used by plume_temporal_evolution() and new_ice_particle_number() as a scaling term.

A seg_ratio value of 1 is the same as not applying any scaling in these two functions.

pycontrails.models.cocip.contrail_properties.temperature_adiabatic_heating(air_temperature, air_pressure, air_pressure_1)

Calculate the ambient air temperature for each waypoint after the wake vortex phase.

This calculation accounts for adiabatic heating.

Parameters:
  • air_temperature (npt.NDArray[np.float64]) – ambient temperature for each waypoint, [\(K\)]

  • air_pressure (npt.NDArray[np.float64]) – initial pressure altitude at each waypoint, before the wake vortex phase, [\(Pa\)]

  • air_pressure_1 (npt.NDArray[np.float64]) – pressure altitude at each waypoint, after the wake vortex phase, [\(Pa\)]

Returns:

npt.NDArray[np.float64] – ambient air temperature after the wake vortex phase, [\(K\)]

Notes

Level 1, see Figure 1 of [Schumann, 2012]

References

pycontrails.models.cocip.contrail_properties.vertical_diffusivity(air_pressure, air_temperature, dT_dz, depth_eff, terminal_fall_speed, sedimentation_impact_factor, eff_heat_rate)

Calculate contrail vertical diffusivity.

Parameters:
  • air_pressure (npt.NDArray[np.float64]) – Pressure altitude at each waypoint, [\(Pa\)]

  • air_temperature (npt.NDArray[np.float64]) – Ambient temperature for each waypoint, [\(K\)]

  • dT_dz (npt.NDArray[np.float64]) – Temperature gradient with respect to altitude (dz), [\(K m^{-1}\)]

  • depth_eff (npt.NDArray[np.float64]) – Effective depth of the contrail plume, [\(m\)]

  • terminal_fall_speed (npt.NDArray[np.float64]) – Terminal fall speed of contrail ice particles, [\(m s^{-1}\)]

  • sedimentation_impact_factor (float) – Enhancement parameter denoted by f_T in eq. (35) Schumann (2012).

  • eff_heat_rate (npt.NDArray[np.float64] | None) – Effective heating rate, i.e., rate of which the contrail plume is heated, [\(K s^{-1}\)]. If None is passed, the radiative heating effects on contrail cirrus properties are not included.

Returns:

npt.NDArray[np.float64] – vertical diffusivity, [\(m^{2} s^{-1}\)]

References

Notes

Accounts for the turbulence-induced diffusive contrail spreading in the vertical direction. See eq. (35) of [Schumann, 2012].

The first term in Eq. (35) of [Schumann, 2012] is (c_V * w’_N^2 / N_BV, where c_V = 0.2 and w’_N^2 = 0.1) is different than outlined below. Here, a constant of 0.01 is used when radiative heating effects are not activated. This update comes from [Schumann and Graf, 2013] , which found that the original formulation estimated thinner contrails relative to satellite observations. The vertical diffusivity was enlarged so that the simulated contrails are more consistent with observations.