as.ctd | R Documentation |
Assemble data into a ctd object. There are two ways this can work.
First, salinity
can be a vector of numeric values, in which case
the other parameters will be interpreted as described below. Second,
salinity
can be an oce object, in which case the action
depends on the object class, as described in the ‘Details’.
as.ctd(
salinity,
temperature = NULL,
pressure = NULL,
conductivity = NULL,
scan = NULL,
time = NULL,
units = NULL,
flags = NULL,
missingValue = NULL,
type = "",
serialNumber = NULL,
ship = NULL,
cruise = NULL,
station = NULL,
startTime = NULL,
longitude = NULL,
latitude = NULL,
deploymentType = "unknown",
pressureAtmospheric = 0,
sampleInterval = NULL,
profile = NULL,
debug = getOption("oceDebug")
)
salinity |
may be (1) a numeric vector holding Practical Salinity,
(2) a list or data frame holding |
temperature |
a numeric vector containing in-situ temperature in
|
pressure |
a numeric vector containing sea pressure values, in decibars.
Typically, this vector has the same length as |
conductivity |
an optional numeric vector containing electrical conductivity ratio through the water column. To convert from raw conductivity in milliSeimens per centimeter divide by 42.914 to get conductivity ratio (see Culkin and Smith, 1980). |
scan |
optional numeric vector holding scan number. If not provided,
this is set to seq_along |
time |
optional vector of times of observation. |
units |
an optional list containing units. If not supplied,
defaults are set for |
flags |
if supplied, this is a list containing data-quality flags. The elements of this list must have names that match the data provided to the object. |
missingValue |
optional missing value, indicating data that should be
taken as |
type |
optional type of CTD, e.g. "SBE" |
serialNumber |
optional serial number of instrument |
ship |
optional string containing the ship from which the observations were made. |
cruise |
optional string containing a cruise identifier. |
station |
optional string containing a station identifier. |
startTime |
optional indication of the start time for the profile,
which is used in some several plotting functions. This is best given as a
POSIXt time, but it may also be a character string
that can be converted to a time with |
longitude |
optional numerical value containing longitude in decimal
degrees, positive in the eastern hemisphere. If this is a single number,
then it is stored in the |
latitude |
similar to |
deploymentType |
character string indicating the type of deployment. Use
|
pressureAtmospheric |
A numerical value (a constant or a vector),
that is subtracted from pressure before storing it in the return value.
(This altered pressure is also used in calculating |
sampleInterval |
optional numerical value indicating the time between samples in the profile. |
profile |
optional positive integer specifying the number of the profile
to extract from an object that has data in matrices, such as for some
|
debug |
an integer specifying whether debugging information is
to be printed during the processing. This is a general parameter that
is used by many |
If the first parameter, salinity
, is an oce object, then
the action depends on the class of that object.
If salinity
is ctd object, then 'as.ctd()1 returns a copy of it.
If salinity
is an argo object, then as.ctd()
calls
argo2ctd()
with that object as its first parameter, along with the value
of profile
and the value of debug
minus 1. All other parameters
provided to as.ctd()
are ignored. Note that Argo notation is retained
in the return value, so that e.g. there is no
metadata item named station
(instead, id
and cycleNumber
are defined),
and no item named startTime
(instead, time
is defined. These name changes
are understood by the summary()
and plot()
functions. Breaking
Change: Until version 1.8-4, as.ctd()
also processed the parameters that are
ignored now. This behaviour was changed because many of those parameters
(e.g. cruise
and ship
) make no sense for Argo data. Users should now
use oceSetMetadata()
to insert additional items as desired.
If salinity
is an rsk object, then as.ctd()
calls rsk2ctd()
with that object as its first argument, along with pressureAtmospheric
,
longitude
, latitude
and debug
minus 1, ignoring all the other
parameters. Note that pressure in the returned object
may need to be adjusted, because rsk
objects may contain either absolute
pressure or sea pressure. This adjustment is handled automatically by
as.ctd
, by examination of the metadata item named pressureType
(described
in the documentation for read.rsk()
). Once the sea pressure is determined,
adjustments may be made with the pressureAtmospheric
argument, although in
that case it is better considered a pressure adjustment than the atmospheric
pressure.
A ctd object.
Dan Kelley, with help from Clark Richards
Culkin, F., and Norman D. Smith, 1980. Determination of the concentration of potassium chloride solution having the same electrical conductivity, at 15 C and infinite frequency, as standard seawater of salinity 35.0000 ppt (Chlorinity 19.37394 ppt). IEEE Journal of Oceanic Engineering, volume 5, pages 22-23.
Other things related to ctd data:
CTD_BCD2014666_008_1_DN.ODF.gz
,
[[,ctd-method
,
[[<-,ctd-method
,
argo2ctd()
,
cnvName2oceName()
,
ctd
,
ctd-class
,
ctd.cnv.gz
,
ctdDecimate()
,
ctdFindProfiles()
,
ctdFindProfilesRBR()
,
ctdRaw
,
ctdRepair()
,
ctdTrim()
,
ctd_aml_type1.csv.gz
,
ctd_aml_type3.csv.gz
,
d200321-001.ctd.gz
,
d201211_0011.cnv.gz
,
handleFlags,ctd-method
,
initialize,ctd-method
,
initializeFlagScheme,ctd-method
,
oceNames2whpNames()
,
oceUnits2whpUnits()
,
plot,ctd-method
,
plotProfile()
,
plotScan()
,
plotTS()
,
read.ctd()
,
read.ctd.aml()
,
read.ctd.itp()
,
read.ctd.odf()
,
read.ctd.odv()
,
read.ctd.saiv()
,
read.ctd.sbe()
,
read.ctd.ssda()
,
read.ctd.woce()
,
read.ctd.woce.other()
,
setFlags,ctd-method
,
subset,ctd-method
,
summary,ctd-method
,
woceNames2oceNames()
,
woceUnit2oceUnit()
,
write.ctd()
library(oce)
# 1. fake data, with default units
pressure <- 1:50
temperature <- 10 - tanh((pressure - 20) / 5) + 0.02 * rnorm(50)
salinity <- 34 + 0.5 * tanh((pressure - 20) / 5) + 0.01 * rnorm(50)
ctd <- as.ctd(salinity, temperature, pressure)
# Add a new column
fluo <- 5 * exp(-pressure / 20)
ctd <- oceSetData(ctd,
name = "fluorescence", value = fluo,
unit = list(unit = expression(mg / m^3), scale = "")
)
summary(ctd)
# 2. fake data, with supplied units (which are the defaults, actually)
ctd <- as.ctd(salinity, temperature, pressure,
units = list(
salinity = list(unit = expression(), scale = "PSS-78"),
temperature = list(unit = expression(degree * C), scale = "ITS-90"),
pressure = list(unit = expression(dbar), scale = "")
)
)
Add the following code to your website.
For more information on customizing the embed code, read Embedding Snippets.