opendrift.models.physics_methods

Attributes

logger

__stokes_coefficients__

__stokes_hs_coefficients__

Classes

PhysicsMethods

Physics methods to be inherited by OpenDriftSimulation class

Functions

wind_drift_factor_from_trajectory(trajectory_dict[, ...])

Estimate wind_drift_fator based on wind and current along given trajectory

distance_between_trajectories(lon1, lat1, lon2, lat2)

Calculate the distances [m] between two trajectories

distance_along_trajectory(lon, lat)

Calculate the distances [m] between points along a trajectory

skillscore_darpa(lon1, lat1, lon2, lat2)

Scoring algorithm used for DARPA float challenge 2021.

skillscore_liu_weissberg(lon_obs, lat_obs, lon_model, ...)

calculate skill score from normalized cumulative seperation

plot_wind_drift_factor(wdf, azimuth[, wmax_plot])

Polar plot of array of wind drift factor, with associated azimuthal offset

oil_wave_entrainment_rate_li2017(dynamic_viscosity, ...)

significant_wave_height_from_wind_neumann_pierson(...)

wave_breaking_fraction_from_wind(wind_speed[, wave_period])

wave_period_from_wind(wind_speed)

verticaldiffusivity_Sundby1983(windspeed, depth[, ...])

Vertical diffusivity from Sundby (1983)

verticaldiffusivity_Large1994(windspeed, depth[, ...])

Vertical diffusivity from Large et al. (1994)

verticaldiffusivity_stepfunction(depth[, MLD, ...])

eddy diffusivity with discontinuity for testing of mixing scheme

gls_tke(windstress, depth, sea_water_density, tke, ...)

From LADIM model.

plot_stokes_profile(profiles[, view, filename])

Plot vertical profile of Stokes drift

stokes_transport_monochromatic(mean_wave_period, ...)

stokes_drift_profile_monochromatic(stokes_u_surface, ...)

Vertical Stokes drift profile assuming a single monochromatic wave

stokes_drift_profile_exponential(stokes_u_surface, ...)

Breivik, O., Janssen, P., Bidlot, J., 2014. Approximate stokes drift profiles in deep water.

stokes_drift_profile_phillips(stokes_u_surface, ...)

Calculate vertical Stokes drift profile from

stokes_drift_profile_windsea_swell(stokes_u_surface, ...)

Calculate vertical Stokes drift profile from

ftle(X, Y, delta, duration)

Calculate Finite Time Lyapunov Exponents

wave_stokes_drift_parameterised(wind, fetch)

Parameterise stokes drift based on pre calculated tables and fetch.

wave_significant_height_parameterised(wind, fetch)

Parameterise significant wave height based on pre calculated tables and fetch.

wind_drag_coefficient(windspeed)

Large and Pond (1981), J. Phys. Oceanog., 11, 324-336.

windspeed_from_stress_polyfit(wind_stress)

Inverting Large and Pond (1981) using polyfit

declination(time)

Solar declination in degrees.

equation_of_time(time)

Equation of time in minutes.

hour_angle(time, longitude)

Solar hour angle in degrees.

solar_elevation(time, longitude, latitude)

Solar elevation in degrees.

Module Contents

opendrift.models.physics_methods.logger
opendrift.models.physics_methods.wind_drift_factor_from_trajectory(trajectory_dict, min_period=None)[source]

Estimate wind_drift_fator based on wind and current along given trajectory

trajectory_dict: dictionary with arrays of same length of the following variables: time, lon, lat, x_wind, y_wind, x_sea_water_velocity, y_sea_water_velocity

Returns array of same length minus one of the fitted wind_drift_factor

opendrift.models.physics_methods.distance_between_trajectories(lon1, lat1, lon2, lat2)[source]

Calculate the distances [m] between two trajectories

opendrift.models.physics_methods.distance_along_trajectory(lon, lat)[source]

Calculate the distances [m] between points along a trajectory assuming the last dimension is the time dimension for the trajectory.

opendrift.models.physics_methods.skillscore_darpa(lon1, lat1, lon2, lat2)[source]

Scoring algorithm used for DARPA float challenge 2021.

Copied from implementation made by Jean Rabault. Assuming 6 positions, separated by 2 days, where first pos (start) are identical.

opendrift.models.physics_methods.skillscore_liu_weissberg(lon_obs, lat_obs, lon_model, lat_model, tolerance_threshold=1)[source]

calculate skill score from normalized cumulative seperation distance. Liu and Weisberg 2011.

Returns the skill score bewteen 0. and 1.

opendrift.models.physics_methods.plot_wind_drift_factor(wdf, azimuth, wmax_plot=None)[source]

Polar plot of array of wind drift factor, with associated azimuthal offset

opendrift.models.physics_methods.oil_wave_entrainment_rate_li2017(dynamic_viscosity, oil_density, interfacial_tension, significant_wave_height=None, wave_breaking_fraction=None, wind_speed=None, sea_water_density=1028.0)[source]
opendrift.models.physics_methods.significant_wave_height_from_wind_neumann_pierson(wind_speed)[source]
opendrift.models.physics_methods.wave_breaking_fraction_from_wind(wind_speed, wave_period=None)[source]
opendrift.models.physics_methods.wave_period_from_wind(wind_speed)[source]
opendrift.models.physics_methods.verticaldiffusivity_Sundby1983(windspeed, depth, mixedlayerdepth=50, background_diffusivity=0)[source]

Vertical diffusivity from Sundby (1983)

  1. Sundby (1983): A one-dimensional model for the vertical

    distribution of pelagic fish eggs in the mixed layer Deep Sea Research (30) pp. 645-661

opendrift.models.physics_methods.verticaldiffusivity_Large1994(windspeed, depth, mixedlayerdepth=50, background_diffusivity=0)[source]

Vertical diffusivity from Large et al. (1994)

Depending on windspeed, depth and mixed layer depth (default 50m).

opendrift.models.physics_methods.verticaldiffusivity_stepfunction(depth, MLD=20, k_above=0.1, k_below=0.02)[source]

eddy diffusivity with discontinuity for testing of mixing scheme

opendrift.models.physics_methods.gls_tke(windstress, depth, sea_water_density, tke, generic_length_scale, gls_parameters=None)[source]

From LADIM model.

opendrift.models.physics_methods.plot_stokes_profile(profiles, view=['vertical', 'birdseye', 'u', 'v'], filename=None)[source]

Plot vertical profile of Stokes drift

Args:
list of dictionary with:

u, v: components of Stokes drift vector z: depth in meters, 0 at surface and negative below (e.g. -5 = 5m depth) kwargs: forwarded to plot method

view:
  • ‘vertical’: magnitude versus depth

  • ‘birdseye’: viewed from above

  • or list of both above with one axis for each case (default)

opendrift.models.physics_methods.stokes_transport_monochromatic(mean_wave_period, significant_wave_height)[source]
opendrift.models.physics_methods.stokes_drift_profile_monochromatic(stokes_u_surface, stokes_v_surface, significant_wave_height, mean_wave_period, z)[source]

Vertical Stokes drift profile assuming a single monochromatic wave Breivik, O., Janssen, P., Bidlot, J., 2014. Approximate stokes drift profiles in deep water. J. Phys. Oceanogr. 44, 2433-2445. doi:10.1175/JPO-D-14-0020.1.

opendrift.models.physics_methods.stokes_drift_profile_exponential(stokes_u_surface, stokes_v_surface, significant_wave_height, mean_wave_period, z)[source]

Breivik, O., Janssen, P., Bidlot, J., 2014. Approximate stokes drift profiles in deep water. J. Phys. Oceanogr. 44, 2433-2445. doi:10.1175/JPO-D-14-0020.1.

opendrift.models.physics_methods.stokes_drift_profile_phillips(stokes_u_surface, stokes_v_surface, significant_wave_height, mean_wave_period, z)[source]

Calculate vertical Stokes drift profile from Breivik et al. 2016, A Stokes drift approximation based on the Phillips spectrum, Ocean Mod. 100

opendrift.models.physics_methods.stokes_drift_profile_windsea_swell(stokes_u_surface, stokes_v_surface, swell_mean_direction_to, swell_mean_period, swell_height, wind_sea_mean_direction_to, wind_sea_mean_period, wind_sea_height, z)[source]

Calculate vertical Stokes drift profile from Breivik, O., and K. H. Christensen, 2020: A Combined Stokes Drift Profile under Swell and Wind Sea.

  1. Phys. Oceanogr., 50, 2819-2833, https://doi.org/10.1175/JPO-D-20-0087.1.

opendrift.models.physics_methods.ftle(X, Y, delta, duration)[source]

Calculate Finite Time Lyapunov Exponents

opendrift.models.physics_methods.__stokes_coefficients__ = None
opendrift.models.physics_methods.wave_stokes_drift_parameterised(wind, fetch)[source]

Parameterise stokes drift based on pre calculated tables and fetch.

opendrift.models.physics_methods.__stokes_hs_coefficients__ = None
opendrift.models.physics_methods.wave_significant_height_parameterised(wind, fetch)[source]

Parameterise significant wave height based on pre calculated tables and fetch.

class opendrift.models.physics_methods.PhysicsMethods[source]

Physics methods to be inherited by OpenDriftSimulation class

static sea_water_density(T=10.0, S=35.0)[source]

The function gives the density of seawater at one atmosphere pressure as given in :

N.P. Fofonoff and R.C. Millard Jr.,1983, Unesco technical papers in marine science no. 44.

S = Salinity in promille of the seawater T = Temperature of the seawater in degrees Celsius

skillscore_trajectory(lon_obs, lat_obs, time_obs, duration=None, max_time_offset=timedelta(seconds=60), method='liu-weissberg', **kwargs)[source]

Calculate skill-score comparing simulated and observed trajectories

A skill score is calculated for each of the modelled trajectories, starting from the closest time/obs of the given modelled trajectory

Parameters

lon_obs : array_like lat_obs : array_like

The longitude and latitudes of the observed trajectory

time_obsarray of datetime

The time corresponding to the observed positions

durationtimedelta

If None, skill score is calculated for the full trajectory, from the start Otherwise, only this sub-part of the trajectory is used.

max_time_offsettimedelta

The maximum allowed time difference between observed and modelled time Default: 1 minute

methodstring

Currently available: ‘liu-wessberg’ To be implemented: ‘darpa’

**kwargs : further arguments for skillscore method (e.g. tolerance_threshold)

Returns

skillscorearray_like

One value for each trajectory

advect_ocean_current(factor=1)[source]
advect_with_sea_ice(factor=1)[source]
advect_wind(factor=1)[source]
stokes_drift(factor=1)[source]
resurface_elements(minimum_depth)[source]
calculate_missing_environment_variables()[source]
wind_speed()[source]
current_speed()[source]
significant_wave_height()[source]
_wave_frequency()[source]
wave_period()[source]
wave_energy()[source]
wave_energy_dissipation()[source]
wave_damping_coefficient()[source]
sea_surface_wave_breaking_fraction()[source]
air_density()[source]
windspeed_from_stress()[source]
solar_elevation()[source]

Solar elevation at present time and position of active elements.

sea_floor_depth()[source]

Sea floor depth (positive) for presently active elements

sea_surface_height()[source]

Sea surface height (positive/negative) for presently active elements

opendrift.models.physics_methods.wind_drag_coefficient(windspeed)[source]

Large and Pond (1981), J. Phys. Oceanog., 11, 324-336.

opendrift.models.physics_methods.windspeed_from_stress_polyfit(wind_stress)[source]

Inverting Large and Pond (1981) using polyfit

opendrift.models.physics_methods.declination(time)[source]

Solar declination in degrees.

opendrift.models.physics_methods.equation_of_time(time)[source]

Equation of time in minutes.

opendrift.models.physics_methods.hour_angle(time, longitude)[source]

Solar hour angle in degrees.

opendrift.models.physics_methods.solar_elevation(time, longitude, latitude)[source]

Solar elevation in degrees.