pycontrails.models.cocip.wake_vortex#

Wave-vortex downwash functions.

Functions

downward_displacement_strongly_stratified(...)

Calculate the maximum contrail downward displacement under strongly stratified conditions.

downward_displacement_weakly_stratified(...)

Calculate the maximum contrail downward displacement under weakly/stably stratified conditions.

effective_time_scale(wingspan, ...)

Calculate the effective time scale of the wake vortex.

initial_contrail_depth(dz_max, ...)

Calculate the initial contrail depth.

initial_contrail_width(wingspan, dz_max)

Calculate the initial contrail width.

max_downward_displacement(wingspan, ...)

Calculate the maximum contrail downward displacement after the wake vortex phase.

normalized_dissipation_rate(epsilon, ...)

Calculate the normalized dissipation rate of the sinking wake vortex.

turbulent_kinetic_energy_dissipation_rate(ds_dz)

Calculate the turbulent kinetic energy dissipation rate (epsilon).

wake_vortex_separation(wingspan)

Calculate the wake vortex separation.

pycontrails.models.cocip.wake_vortex.downward_displacement_strongly_stratified(wingspan, true_airspeed, aircraft_mass, rho_air, n_bv)#

Calculate the maximum contrail downward displacement under strongly stratified conditions.

Parameters:
  • wingspan (npt.NDArray[np.float_] | float) – aircraft wingspan, [\(m\)]

  • true_airspeed (npt.NDArray[np.float_]) – true airspeed for each waypoint, [\(m s^{-1}\)]

  • aircraft_mass (npt.NDArray[np.float_] | float) – aircraft mass for each waypoint, [\(kg\)]

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

  • n_bv (npt.NDArray[np.float_]) – Brunt-Vaisaila frequency, [\(s^{-1}\)]

Returns:

npt.NDArray[np.float_] – Maximum contrail downward displacement, strongly stratified conditions, [\(m\)]

Notes

See section 2.5 (pg 547 - 548) of [Schumann, 2012].

References

pycontrails.models.cocip.wake_vortex.downward_displacement_weakly_stratified(wingspan, true_airspeed, aircraft_mass, rho_air, n_bv, dz_max_strong, ds_dz, t_0, effective_vertical_resolution, wind_shear_enhancement_exponent)#

Calculate the maximum contrail downward displacement under weakly/stably stratified conditions.

Parameters:
  • wingspan (npt.NDArray[np.float_] | float) – aircraft wingspan, [\(m\)]

  • true_airspeed (npt.NDArray[np.float_]) – true airspeed for each waypoint, [\(m s^{-1}\)]

  • aircraft_mass (npt.NDArray[np.float_] | float) – aircraft mass for each waypoint, [\(kg\)]

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

  • n_bv (npt.NDArray[np.float_]) – Brunt-Vaisaila frequency, [\(s^{-1}\)]

  • dz_max_strong (npt.NDArray[np.float_]) – Max contrail downward displacement under strongly stratified conditions, [\(m\)]

  • ds_dz (npt.NDArray[np.float_]) – Difference in wind speed over dz in the atmosphere, [\(m s^{-1} / m\)]

  • t_0 (npt.NDArray[np.float_]) – Wake vortex effective time scale, [\(s\)]

  • effective_vertical_resolution (float) – Passed through to wind_shear.wind_shear_enhancement_factor(), [\(m\)]

  • wind_shear_enhancement_exponent (npt.NDArray[np.float_] | float) – Passed through to wind_shear.wind_shear_enhancement_factor()

Returns:

npt.NDArray[np.float_] – Maximum contrail downward displacement, weakly/stably stratified conditions, [\(m\)]

Notes

See section 2.5 (pg 548) of [Schumann, 2012].

References

pycontrails.models.cocip.wake_vortex.effective_time_scale(wingspan, true_airspeed, aircraft_mass, rho_air)#

Calculate the effective time scale of the wake vortex.

Parameters:
  • wingspan (npt.NDArray[np.float_]) – aircraft wingspan, [\(m\)]

  • true_airspeed (npt.NDArray[np.float_]) – true airspeed for each waypoint, [\(m \ s^{-1}\)]

  • aircraft_mass (npt.NDArray[np.float_]) – aircraft mass for each waypoint, [\(kg\)]

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

Returns:

npt.NDArray[np.float_] – Wake vortex effective time scale, [\(s\)]

Notes

See section 2.5 (pg 547) of [Schumann, 2012].

References

pycontrails.models.cocip.wake_vortex.initial_contrail_depth(dz_max, initial_wake_vortex_depth)#

Calculate the initial contrail depth.

Parameters:
  • dz_max (npt.NDArray[np.float_]) – Max contrail downward displacement after the wake vortex phase, [\(m\)]

  • initial_wake_vortex_depth (float | npt.NDArray[np.float_]) – Initial wake vortex depth scaling factor. Denoted C_D0 in eq (14) in [Schumann, 2012].

Returns:

npt.NDArray[np.float_] – Initial contrail depth, [\(m\)]

pycontrails.models.cocip.wake_vortex.initial_contrail_width(wingspan, dz_max)#

Calculate the initial contrail width.

Parameters:
  • wingspan (npt.NDArray[np.float_] | float) – aircraft wingspan, [\(m\)]

  • dz_max (npt.NDArray[np.float_]) – Max contrail downward displacement after the wake vortex phase, [\(m\)] Only the size of this array is used; the values are ignored.

Returns:

npt.NDArray[np.float_] – Initial contrail width, [\(m\)]

pycontrails.models.cocip.wake_vortex.max_downward_displacement(wingspan, true_airspeed, aircraft_mass, air_temperature, dT_dz, ds_dz, air_pressure, effective_vertical_resolution, wind_shear_enhancement_exponent)#

Calculate the maximum contrail downward displacement after the wake vortex phase.

Parameters:
  • wingspan (npt.NDArray[np.float_] | float) – aircraft wingspan, [\(m\)]

  • true_airspeed (npt.NDArray[np.float_]) – true airspeed for each waypoint, [\(m s^{-1}\)]

  • aircraft_mass (npt.NDArray[np.float_] | float) – aircraft mass for each waypoint, [\(kg\)]

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

  • dT_dz (npt.NDArray[np.float_]) – potential temperature gradient, [\(K m^{-1}\)]

  • ds_dz (npt.NDArray[np.float_]) – Difference in wind speed over dz in the atmosphere, [\(m s^{-1} / m\)]

  • air_pressure (npt.NDArray[np.float_]) – pressure altitude at each waypoint, [\(Pa\)]

  • effective_vertical_resolution (float) – Passed through to wind_shear.wind_shear_enhancement_factor(), [\(m\)]

  • wind_shear_enhancement_exponent (npt.NDArray[np.float_] | float) – Passed through to wind_shear.wind_shear_enhancement_factor()

Returns:

npt.NDArray[np.float_] – Max contrail downward displacement after the wake vortex phase, [\(m\)]

References

pycontrails.models.cocip.wake_vortex.normalized_dissipation_rate(epsilon, wingspan, true_airspeed, aircraft_mass, rho_air)#

Calculate the normalized dissipation rate of the sinking wake vortex.

Parameters:
  • epsilon (npt.NDArray[np.float_]) – turbulent kinetic energy dissipation rate, [\(m^{2} s^{-3}\)]

  • wingspan (npt.NDArray[np.float_] | float) – aircraft wingspan, [\(m\)]

  • true_airspeed (npt.NDArray[np.float_]) – true airspeed for each waypoint, [\(m s^{-1}\)]

  • aircraft_mass (npt.NDArray[np.float_]) – aircraft mass for each waypoint, [\(kg\)]

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

Returns:

npt.NDArray[np.float_] – Normalized dissipation rate of the sinking wake vortex

Notes

See page 548 of [Schumann, 2012].

References

pycontrails.models.cocip.wake_vortex.turbulent_kinetic_energy_dissipation_rate(ds_dz, shear_enhancement_factor=1.0)#

Calculate the turbulent kinetic energy dissipation rate (epsilon).

The shear enhancement factor is used to account for any sub-grid scale turbulence.

Parameters:
  • ds_dz (npt.NDArray[np.float_]) – Difference in wind speed over dz in the atmosphere, [\(m s^{-1} / m\)]

  • shear_enhancement_factor (npt.NDArray[np.float_] | float) – Multiplication factor to enhance the wind shear

Returns:

npt.NDArray[np.float_] – turbulent kinetic energy dissipation rate, [\(m^{2} s^{-3}\)]

Notes

  • See eq. (37) in [Schumann, 2012].

  • In a personal correspondence, Dr. Schumann identified a print error in Eq. (37) of the 2012 paper where the shear term should not be squared. The correct equation is listed in Eq. (13) [Schumann and Gerz, 1995].

References

pycontrails.models.cocip.wake_vortex.wake_vortex_separation(wingspan)#

Calculate the wake vortex separation.

Parameters:

wingspan (npt.NDArray[np.float_] | float) – aircraft wingspan, [\(m\)]

Returns:

npt.NDArray[np.float_] – wake vortex separation, [\(m\)]