PhysOcean
Sea-water properties
PhysOcean.density
— Methoddensity(S,T,p)
Compute the density of sea-water (kg/m³) at the salinity S
(psu, PSS-78), temperature T
(degree Celsius, ITS-90) and pressure p
(decibar) using the UNESCO 1983 polynomial.
Fofonoff, N.P.; Millard, R.C. (1983). Algorithms for computation of fundamental properties of seawater. UNESCO Technical Papers in Marine Science, No. 44. UNESCO: Paris. 53 pp. http://web.archive.org/web/20170103000527/http://unesdoc.unesco.org/images/0005/000598/059832eb.pdf
PhysOcean.secant_bulk_modulus
— Methodsecant_bulk_modulus(S,T,p)
Compute the secant bulk modulus of sea-water (bars) at the salinity S
(psu, PSS-78), temperature T
(degree Celsius, ITS-90) and pressure p
(decibar) using the UNESCO polynomial 1983.
Fofonoff, N.P.; Millard, R.C. (1983). Algorithms for computation of fundamental properties of seawater. UNESCO Technical Papers in Marine Science, No. 44. UNESCO: Paris. 53 pp. http://web.archive.org/web/20170103000527/http://unesdoc.unesco.org/images/0005/000598/059832eb.pdf
PhysOcean.freezing_temperature
— Methodfreezing_temperature(S)
Compute the freezing temperature (in degree Celsius) of sea-water based on the salinity S
(psu).
PhysOcean.potential_temperature
— Functionpotential_temperature(S,T,P,PR)
Potential temperature (degree Celsius, ITS-90) as per UNESCO 1983 report of a water mass with the salinity S
(psu, PSS-78) and temperature T
(degree Celsius, ITS-90) at the pressure P
(db) relative to the reference pressure PR
(db).
PhysOcean.potential_density
— Functionpotential_density(S,T,P,PR)
Potential density (kg/m^3) of a water mass with the salinity S
(psu, PSS-78) and temperature T
(degree Celsius, ITS-90) at the pressure P
(db) relative to the reference pressure PR
(db).
PhysOcean.N²
— FunctionPhysOcean.N²(S,T,P,latitude)
Brunt-Väisälä or buoyancy frequency squared (1/s²) at mid-depth for a profile with the salinity S
(psu, PSS-78), temperature T
(degree Celsius, ITS-90) and pressure P
(decibar). The latitude lat
(degrees) is used to compute the acceleration due to gravity.
\[N^2 = - \frac{g}{ρ} \frac{∂ρ}{∂z}\]
where z is negative in water.
Heat fluxes
PhysOcean.latentflux
— Methodlatentflux(Ts,Ta,r,w,pa)
Compute the latent heat flux (W/m²) using the sea surface temperature Ts
(degree Celsius), the air temperature Ta
(degree Celsius), the relative humidity r
(0 ≤ r ≤ 1, pressure ratio, not percentage), the wind speed w
(m/s) and the air pressure (hPa).
PhysOcean.longwaveflux
— Methodlongwaveflux(Ts,Ta,e,tcc)
Compute the long wave heat flux (W/m²) using the sea surface temperature Ts
(degree Celsius), the air temperature Ta
(degree Celsius), the wate vapour pressure e
(hPa) and the total cloud coverage ttc
(0 ≤ tcc ≤ 1).
PhysOcean.sensibleflux
— Methodsensibleflux(Ts,Ta,w)
Compute the sensible heat flux (W/m²) using the wind speed w
(m/s), the sea surface temperature Ts
(degree Celsius), the air temperature Ta
(degree Celsius).
PhysOcean.solarflux
— Methodsolarflux(Q,al)
Compute the solar heat flux (W/m²)
PhysOcean.vaporpressure
— Methodvaporpressure(T)
Compute vapour pressure of water at the temperature T
(degree Celsius) in hPa using Tetens equations. The temperature must be postive.
Monteith, J.L., and Unsworth, M.H. 2008. Principles of Environmental Physics. Third Ed. AP, Amsterdam. http://store.elsevier.com/Principles-of-Environmental-Physics/John-Monteith/isbn-9780080924793
Filtering
PhysOcean.gausswin
— Functiongausswin(N, α = 2.5)
Return a Gaussian window with N
points with a standard deviation of (N-1)/(2 α).
PhysOcean.gaussfilter
— Methodgaussfilter(x,N)
Filter the vector x
with a N
-point Gaussian window.
Earth
PhysOcean.coriolisfrequency
— Methodcoriolisfrequency(latitude)
Provides coriolisfrequency et given latidudes in DEGREES from -90 Southpole to +90 Northpole
PhysOcean.earthgravity
— Methodearthgravity(latitude)
Provides gravity in m/s2 at ocean surface at given latidudes in DEGREES from -90 Southpole to +90 Northpole
GFD
PhysOcean.integraterhoprime
— Methodrhoi = integraterhoprime(rhop,z,dim=ndims(rhop))
Integrates density anomalies over depth. When used with gravity, assuming gravity is independant on z, it can be used to calculate dynamic pressure up to a constant. Function can be used with 1D, 2D, ...
Input:
rhop
: density anomaly arrayz
: vertical position array. Zero at surface and positive downward, same dimensions as rhopdim
: along which dimension depth is found and integral is performed. If not provided, last dimension is taken
Output:
rhoi
: Integrated value to the same levels as on which rhop where given. So basically total density anomaly ABOVE the current depth
Note:
Compute vertical integral of density anomalies
PhysOcean.stericheight
— Functionssh=stericheight(rhoi,z,zlevel,dim::Integer=ndims(rhoi))
Input:
rhoi
: integrated density anomalies (from a call to integraterhoprime)z
: array of vertical positionszlevel
: integer for the zlevel on which no motion is assumeddim
: along which dimension depth is found . If not provided last dimension is used
Output:
ssh
: steric height. space dimensions as for rhoi in which direction dim is taken out
Compute steric height with respect to given depth level presently provided as index , not depth
PhysOcean.geostrophy
— Methodvelocity,ssh,fluxes=geostrophy(mask::BitArray,rhop,pmnin,xiin;dim::Integer=ndims(rhop),ssh=(),znomotion=0,fillin=true)
Input:
mask
: Boolean array with true in water and false on land.rhop
: density anomaly (rho-1025) array on the same grid as the mask.pmnin
: tuple of metrics as in divand, but to get velocities in m/s the metrics need to be in per meters too.xiin
: tuple position of the grid points.dim
: optional paramter telling which index in the arrays corresponds to the vertical direction. By default the last dimension is used.ssh
: array as mask for which the vertical direction is taken out. Corresponds to sea surface height in meters. Default is no used but diagnosedznomotion
: index in the vertical direction where a level of no motion is assumedfillin
: Boolean telling if a filling of land points using water points at the same level is to be used. Default is yes.
Output:
velocity
tuple of velocity components NORMAL and to the left of each coordinate lineeta
: sea surface height deduced. If a ssh was provided in input it returs ssh but filled in on land.fluxes
: integrated velocities across sections in each horizontal direction. Same conventions as for velocities
Note:
Calculates geostrophic velocities. Works with one or two horizontal dimensions and additional (time) dimensions. Dimensions are supposed to be ordered horizontal, vertical, other dimensions You must either provide the index for the level of no motion or ssh eta. NOTE THAT THE LEVEL IS AN INDEX NUMBER FOR THE MOMENT Dimensions of ssh must be the same as rhop in which vertical dimension has been taken out If you force fillin=false, then you must have created the density array without missing values outside of this call, as well as ssh if you provide it.
PhysOcean.streamfunctionvolumeflux
— Methodpsifluxes=streamfunctionvolumeflux(mask::BitArray,velocities,pmnin,xiin;dim::Integer=ndims(mask))
Input:
mask
: Boolean array with true in water and false on land.velocities
: tuple of arrays on the same grid as the mask. Each tuple element is a velocity field normal to the corresponding direction in spacepmnin
: tuple of metrics as in divand, but to get velocities in m/s the metrics need to be in per meters too.xiin
: tuple position of the grid points.dim
: optional paramter telling which index in the arrays corresponds to the vertical direction. By default the last dimension is used.
Output:
psifluxes
tuple of volume fluxes at each depth and direction NORMAL and to the left of each coordinate line
Calculates volume flux streamfunction calculated from the surface. The value of this field provides the total flow (in Sverdrup) across the section above the depth of the zlevel looked at.
Tides
PhysOcean.ap2ep
— Functionsemimajor, eccentricity, inclination, phase = ap2ep(u::Complex,v::Complex)
semimajor, eccentricity, inclination, phase = ap2ep(amplitude_u,phase_u,amplitude_v,phase_v)
Return the tidal elipise parameter semi-major axis (semimajor
), eccentricity
, inclination (degrees, angle between east and the maximum current) and phase (degrees) from the amplitude and phase of the meridional and zonal velocity. The velocities are related to the amplitude and phase by:
u(t) = aᵤ cos(ω t - ϕᵤ) v(t) = aᵥ cos(ω t - ϕᵥ)
where ω is the tidal angular frequency and t the time.
The code is based on the technical report Ellipse Parameters Conversion and Velocity Profile For Tidal Currents from Zhigang Xu.
semiminor = semimajor * eccentricity
Example
using PyPlot, PhysOcean
amplitude_u = 1.3
amplitude_v = 1.12
phase_u = 190
phase_v = 350
semimajor, eccentricity, inclination, phase = ap2ep(amplitude_u,phase_u,amplitude_v,phase_v)
semiminor = semimajor * eccentricity
phase_ = 0:1:360
u = amplitude_u * cosd.(phase_ .- phase_u)
v = amplitude_v * cosd.(phase_ .- phase_v)
plot(u,v,"-",label="tidal ellipse")
plot([0,semimajor * cosd(inclination)],[0,semimajor * sind(inclination)],"-",label = "semi-major axis")
plot([0,semiminor * cosd(inclination+90)],[0,semiminor * sind(inclination+90)],"-",label = "semi-minor axis")
legend()
PhysOcean.ep2ap
— Functionamplitude_u,phase_u,amplitude_v,phase_v = ep2ap(semimajor, eccentricity, inclination, phase)
Inverse of ap2eps
.
PhysOcean.rayleigh_criterion
— FunctionTmin = PhysOcean.rayleigh_criterion(f1,f2)
Compute the minimum duration Tmin
to separate the frequency f1
and f2
following the Rayleigh criterion. Tmin
has the inverse of the units of f1
and f2
.
For example to resolve the M2 and S2 tides, the duration of the time series has to be at least 14.765 days:
using PhysOcean
f_M2 = 12.4206 # h⁻¹
f_S2 = 12 # h⁻¹
Tmin = PhysOcean.rayleigh_criterion(f_M2,f_S2)
Tmin/24
# output 14.765
Reference:
Bruce B. Parker, Tidal Analysis and Prediction, NOAA Special Publication NOS CO-OPS 3, page 84
Data download
CMEMS
PhysOcean.CMEMS.download
— FunctionCMEMS.download(lonr,latr,timerange,param,username,password,basedir[; indexURLs = ...])
The data service has been discontinued by CMEMS and replaced by a python-specific API.
Download in situ data within the longitude range lonr
(an array or tuple with two elements: the minimum longitude and the maximum longitude), the latitude range latr
(likewise), time range timerange
(an array or tuple with two DateTime
structures: the starting date and the end date) from the CMEMS (Copernicus Marine environment monitoring service) in situ service [1]. param
is one of the parameter codes as defined in [2] or [3]. username
and password
are the credentials to access data [1] and basedir
is the directory under which the data is saved. indexURLs
is a list of the URL to the index_history.txt
file. Per default, it includes the URLs of the Baltic, Arctic, North West Shelf, Iberian, Mediteranean and Black Sea Thematic Assembly Center.
As these URLs might change, the latest version of the URLs to the indexes can be obtained at [1].
Example
julia> username = "..."
julia> password = "..."
julia> lonr = [7.6, 12.2]
julia> latr = [42, 44.5]
julia> timerange = [DateTime(2016,5,1),DateTime(2016,8,1)]
julia> param = "TEMP"
julia> basedir = "/tmp"
julia> files = CMEMS.download(lonr,latr,timerange,param,username,password,basedir)
PhysOcean.CMEMS.load
— Functiondata,lon,lat,z,time,ids = load(T,fname::TS,param; qualityflags = [good_data, probably_good_data]) where TS <: AbstractString
obsdata,obslon,obslat,obsz,obstime,obsids = CMEMS.load(T,fnames,param;
qualityflags = ...,; prefixid = "")
Load all data in the vector of file names fnames
corresponding to the parameter param
as the data type T
. Only the data with the quality flags CMEMS.good_data
and CMEMS.probably_good_data
are loaded per default. The output parameters correspondata to the data, longitude, latitude, depth, time (as DateTime
) and an identifier (as String
).
If prefixid
is specified, then the observations identifier are prefixed with prefixid
.
See also CMEMS.download
.
World Ocean Database (US NODC)
PhysOcean.WorldOceanDatabase.download
— Functiondirnames,indexnames = WorldOceanDatabase.download(lonrange,latrange,timerange,
variable,email,basedir)
Download data using the NODC web-service. The range parameters are vectors from with the frist element is the lower bound and the last element is the upper bound. The parameters of the functions will be transmitted to nodc.noaa.gov (http://www.noaa.gov/privacy.html). Note that no XBT corrections are applied. The table below show the avialable variable and their units.
Example:
dirnames,indexnames = WorldOceanDatabase.download([0,10],[30,40], [DateTime(2000,1,1),DateTime(2000,2,1)], "Temperature","your@email.com,"/tmp")
Variables | Unit |
---|---|
Temperature | °C |
Salinity | unitless |
Oxygen | ml l⁻¹ |
Phosphate | µM |
Silicate | µM |
Nitrate and Nitrate+Nitrite | µM |
pH | unitless |
Chlorophyll | µg l⁻¹ |
Plankton | multiple |
Alkalinity | meq l⁻¹ |
Partial Pressure of Carbon Dioxide | µatm |
Dissolved Inorganic Carbon | mM |
Transmissivity | m⁻¹ |
Pressure | dbar |
Air temperature | °C |
CO2 warming | °C |
CO2 atmosphere | ppm |
Air pressure | mbar |
Tritium | TU |
Helium | nM |
Delta Helium-3 | % |
Delta Carbon-14 | ᵒ/ᵒᵒ |
Delta Carbon-13 | ᵒ/ᵒᵒ |
Argon | nM |
Neon | nM |
Chlorofluorocarbon 11 (CFC 11) | pM |
Chlorofluorocarbon 12 (CFC 12) | pM |
Chlorofluorocarbon 113 (CFC 113) | pM |
Delta Oxygen-18 | ᵒ/ᵒᵒ |
PhysOcean.WorldOceanDatabase.load
— Function load(T,dirname,indexname,varname)
Load all profiles with the NetCDF variable varname
in dirname
indexed with the NetCDF file indexname
. T is the type (e.g. Float64) for numeric return values.
value,lon,lat,depth,obstime,id = load(T,dirnames,varname)
Load a list of directories dirnames
.
obsvalue,obslon,obslat,obsdepth,obstime,obsid = WorldOceanDatabase.load(T,basedir::AbstractString,varname; prefixid = "")
Load a list profiles under the directory basedir
assuming basedir
was populated by WorldOceanDatabase.download
. If the prefixid
is specified, then the observations identifier are prefixed with prefixid
.
Example:
basedir = expanduser("~/Downloads/woddownload")
varname = "Temperature"
prefixid = "1977-"
obsvalue,obslon,obslat,obsdepth,obstime,obsid = WorldOceanDatabase.load(Float64,basedir,varname; prefixid = prefixid);
ARGO
PhysOcean.ARGO.load
— Functionobsdata,obslon,obslat,obsz,obstime,obsids = ARGO.load(T,fnames,param;
qualityflags = ...,; prefixid = "")
Load all data in the vector of file names fnames
corresponding to the parameter param
as the data type T
. Only the data with the quality flags ARGO.good_data
and ARGO.probably_good_data
are loaded per default. The output parameters correspondata to the data, longitude, latitude, depth, time (as DateTime
) and an identifier (as String
).
If prefixid
is specified, then the observations identifier are prefixed with prefixid
.
Example
using PhysOcean, Glob
# directory containing only ARGO/CORA netCDF files, e.g.
# <basedir>/CORA-5.2-mediterrane-1950/CO_DMQCGL01_19500103_PR_SH.nc
basedir = "/some/dir"
filenames = glob("*/*nc",basedir)
# load the data as double precision
T = Float64
# name of the variable in the NetCDF files
varname = "TEMP"
# EDMO code of Coriolis
prefixid = "4630-"
obsvalue,obslon,obslat,obsdepth,obstime,obsids = ARGO.load(T,filenames,varname;
prefixid = prefixid)
See also CMEMS.load
.