swDynamicHeight {oce} | R Documentation |
Dynamic Height of a Seawater Profile
Description
Compute the dynamic height of a column of seawater.
Usage
swDynamicHeight(
x,
referencePressure = 2000,
subdivisions = 500,
rel.tol = .Machine$double.eps^0.25,
eos = getOption("oceEOS", default = "gsw")
)
Arguments
x |
a section object. |
referencePressure |
reference pressure (dbar). If this exceeds the
highest pressure supplied to |
subdivisions |
number of subdivisions for call to
|
rel.tol |
absolute tolerance for call to |
eos |
equation of state, either |
Details
If the first argument is a section
, then dynamic height is calculated
for each station within a section, and returns a list containing distance
along the section along with dynamic height.
If the first argument is a ctd
, then this returns just a single
value, the dynamic height.
If eos="unesco"
, processing is as follows. First, a piecewise-linear
model of the density variation with pressure is developed using
stats::approxfun()
. (The option rule=2
is used to
extrapolate the uppermost density up to the surface, preventing a possible a
bias for bottle data, in which the first depth may be a few metres below the
surface.) A second function is constructed as the density of water with
salinity 35PSU, temperature of 0C, and pressure as in the
ctd
. The difference of the reciprocals of these densities, is then
integrated with stats::integrate()
with pressure limits 0
to referencePressure
. (For improved numerical results, the variables
are scaled before the integration, making both independent and dependent
variables be of order one.)
If eos="gsw"
, gsw::gsw_geo_strf_dyn_height()
is used
to calculate a result in m^2/s^2, and this is divided by
9.7963.
If pressures are out of order, the data are sorted. If any pressure
is repeated, only the first level is used.
If there are under 4 remaining distinct
pressures,
NA
is returned, with a warning.
Value
In the first form, a list containing distance
, the distance
(km( from the first station in the section and height
, the dynamic
height (m). In the second form, a single value, containing the
dynamic height (m).
Sample of Usage
library(oce) data(section) # Dynamic height and geostrophy par(mfcol=c(2, 2)) par(mar=c(4.5, 4.5, 2, 1)) # Left-hand column: whole section # (The smoothing lowers Gulf Stream speed greatly) westToEast <- subset(section, 1<=stationId&stationId<=123) dh <- swDynamicHeight(westToEast) plot(dh$distance, dh$height, type="p", xlab="", ylab="dyn. height [m]") ok <- !is.na(dh$height) smu <- supsmu(dh$distance, dh$height) lines(smu, col="blue") f <- coriolis(section[["station", 1]][["latitude"]]) g <- gravity(section[["station", 1]][["latitude"]]) v <- diff(smu$y)/diff(smu$x) * g / f / 1e3 # 1e3 converts to m plot(smu$x[-1], v, type="l", col="blue", xlab="distance [km]", ylab="velocity (m/s)") # right-hand column: gulf stream region, unsmoothed gs <- subset(section, 102<=stationId&stationId<=124) dh.gs <- swDynamicHeight(gs) plot(dh.gs$distance, dh.gs$height, type="b", xlab="", ylab="dyn. height [m]") v <- diff(dh.gs$height)/diff(dh.gs$distance) * g / f / 1e3 plot(dh.gs$distance[-1], v, type="l", col="blue", xlab="distance [km]", ylab="velocity (m/s)")
Author(s)
Dan Kelley
References
Gill, A.E., 1982. Atmosphere-ocean Dynamics, Academic Press, New York, 662 pp.
See Also
Other functions that calculate seawater properties:
T68fromT90()
,
T90fromT48()
,
T90fromT68()
,
computableWaterProperties()
,
locationForGsw()
,
swAbsoluteSalinity()
,
swAlphaOverBeta()
,
swAlpha()
,
swBeta()
,
swCSTp()
,
swConservativeTemperature()
,
swDepth()
,
swLapseRate()
,
swN2()
,
swPressure()
,
swRho()
,
swRrho()
,
swSCTp()
,
swSR()
,
swSTrho()
,
swSigma0()
,
swSigma1()
,
swSigma2()
,
swSigma3()
,
swSigma4()
,
swSigmaTheta()
,
swSigmaT()
,
swSigma()
,
swSoundAbsorption()
,
swSoundSpeed()
,
swSpecificHeat()
,
swSpice()
,
swSstar()
,
swTFreeze()
,
swTSrho()
,
swThermalConductivity()
,
swTheta()
,
swViscosity()
,
swZ()