OceanAnalysis Functions

Data types:

OceanAnalysis.OAType
Base type in the OceanAnalysis package.

This is an abstract type, from which other types in the library are derived.

source
OceanAnalysis.AdpType
A type to hold acoustic-Doppler profiler data

Adp is a type used to store data from an ADP (acoustic-Doppler profiler). It holds two items: data (a Dict that holds data in array form) and metadata (a Dict with information about the data).

Adp objects may be created with read_adp_rdi, which only handles RDI data. (Use the oce R package if you need other instrument types.)

See the online documentation for an illustration.

source
OceanAnalysis.AmsrType
A type to hold AMSR satellite data

This holds AMSR satellite data as read by read_amsr.

The metadata element is a Dict that holds the source filename, the sensor, the name of the stored variable, the observation interval ("time_coverage_start" and time_coverage_end) and the vectors (longitude and latitude) that define the grid. The data element is a matrix holding the gridded data.

See the online documentation for an illustration.

source
OceanAnalysis.ArgoType
A type to hold Argo data

Argo is a type used to store data from Argo floats.

Argo objects hold two components: a DataFrame named data that holds the actual data, and a Dict named metadata that stores information about the data, e.g. the sampling time, the location, etc. Objects of type Argo are returned by read_argo, which reads the NetCDF files in which Argo data are distributed.

For many purposes, it may be useful to convert Argo objects into Ctd objects, using as_ctd.

See the online documentation for an illustration.

source
OceanAnalysis.CtdType
A type to hold CTD data

Ctd is a type used to store data from CTD instruments.

Ctd objects hold two items: (1) a DataFrame named data that holds the actual data, including pressure, salinity, and temperature, perhaps along with other data and (2) a Dict named metadata that stores information about the data, such as the time of observation and the latitude and longitude at which the observation was made.

Objects of type Ctd are returned by as_ctd, read_ctd_cnv, read_ctd_rsk and read_argo. Such objects can be passed to plotting functions plot_profile and plot_TS, and to some functions relating to seawater properties, such as SA and other TEOS-10 related functions, as well as functions relating to the distributions of such properties, such as N2.

See the online documentation for an illustration.

source
OceanAnalysis.CoastlineType
A type to hold coastline data

Coastline is a type to hold coastline data. Its metadata element is a Dict that may hold the source filename or other information (or may be empty). Its data element is a DataFrame holding columns named longitude and latitude, with NaN values separating continents and/or elements within them such as countries or subregions of countries.

See the online documentation for an illustration.

source
OceanAnalysis.NonnaType
A type to hold NONNA bathymetry data

This holds bathymetric data as read by read_nonna.

The metadata element is a Dict that holds the source filename along with vectors named longitude and latitude that define the grid. The data element stores a matrix of bathymetry data in terms of height above mean sea level, in metres.

See the online documentation for an illustration.

source
OceanAnalysis.TopographyType
A type to hold topography data

This holds topography data as read by read_topography.

The metadata element is a Dict that holds the source filename along with vectors holding the longitude and latitude values that define the grid.

The data element stores a matrix of topography data in terms of height above mean sea level, in metres.

See the online documentation for an illustration.

source

Functions:

OceanAnalysis.argo_id_cycleFunction
Split Argo "id_cycle" into components id and cycle

Examples

julia> using OceanAnalysis

julia> argo_id_cycle("4902911_095")
2-element Vector{SubString{String}}:
 "4902911"
 "095"
source
OceanAnalysis.as_ctdFunction
as_ctd(a::Argo; add_teos::Bool=false, debug::Int64=0)

Convert an Argo object into a Ctd object.

Return Value

This returns a Ctd object, with metadata and data copied from a, and possibly with new data columns holding computed values of some key TEOS-10 values.

Arguments

  • a an [Argo] object.

Keywords

  • add_teos a logical value indicating whether to add TEOS-10 items (e.g. SA) to the data portion of the return value.

  • debug: an optional value that, if it exceeds 0, indicates that debugging output should be printed during processing.

source
as_ctd(salinity::Union{AbstractVector,AbstractRange},
    temperature::Union{AbstractVector,AbstractRange},
    pressure::Union{AbstractVector,AbstractRange};
    longitude::Real=-63.0, latitude::Real=45.0, time=nothing,
    add_teos::Bool=true, debug::Int64=0)

Construct a Ctd object, given S, T, p, and possibly a location.

Returns a Ctd object with a data element that is a data frame holding the provided water properties, along with computed Absolute Salinity (SA) Conservative Temperature (CT), potential density anomaly relative to the surface pressure (sigma0) and spiciness with respect to surface pressure (spiciness0). The object also holds a metadata element that holds longitude, latitude and time. If either longitude or latitude is NaN, thenSA, etc. are computed assuming a mid-Atlantic location (-30E and 30N).

Arguments

  • salinity: measured salinity values, in Practical Salinity units.

  • temperature: measured temperature values, in degrees Celsius.

  • pressure: measured sea pressure values, in dbar.

Keywords

  • longitude: observation longitude, in degrees East. The default is a location in the North Atlantic).

  • latitude: observation latitude, in degrees North. The default is a location in the North Atlantic).

  • time: an optional indication of the measurement start time.

  • add_teos: an optional indication of whether to add SA, CT, sigma0 and spiciness0 to the data component of the return value.

  • debug: an optional value that, if it exceeds 0, indicates that debugging output should be printed during processing.

Examples

using OceanAnalysis
julia> using OceanAnalysis

julia> as_ctd([32.], [15.], [0.], add_teos=false)
Ctd(Dict{String, Any}("filename" => nothing, "latitude" => [45.0], "time" => nothing, "longitude" => [-63.0]), 1×3 DataFrame
 Row │ salinity  temperature  pressure
     │ Float64   Float64      Float64
─────┼─────────────────────────────────
   1 │     32.0         15.0       0.0)

julia> as_ctd([32.], [15.], [0.], longitude=40.0, latitude=-63.0)
Ctd(Dict{String, Any}("filename" => nothing, "latitude" => [-63.0], "time" => nothing, "longitude" => [40.0]), 1×7 DataFrame
 Row │ salinity  temperature  pressure  SA       CT       sigma0   spiciness0
     │ Float64   Float64      Float64   Float64  Float64  Float64  Float64
─────┼────────────────────────────────────────────────────────────────────────
   1 │     32.0         15.0       0.0  32.1539  15.0641  23.6671   0.0704222)
source
OceanAnalysis.as_sectionFunction
as_section(ctds::Vector{Ctd}; name::String="", source::String="", debug::Int64=0)

Construct a Section object, which olds one or more Ctd objects.

Returns a Section object with a data element that is a vector of Ctd objects, along with a metadata element that holds the name of the section and the source of the data.

Arguments

  • ctds: Vector of Ctd objects.

Keywords

  • name: the name of the section, often an 'EXPOCODE'.

  • source: an indication of the data source, often a URL.

  • debug: an optional value that, if it exceeds 0, indicates that debugging output should be printed during processing.

Examples

julia> using OceanAnalysis
julia> # Create fake data
julia> f = joinpath(dirname(dirname(pathof(OceanAnalysis))), "data", "ctd.cnv");
julia> a = read_ctd_cnv(f, add_teos=false);
julia> b = read_ctd_cnv(f, add_teos=false);
julia> b.data.salinity = 1.0 .+ b.data.salinity;
julia> s = as_section([a; b], name="Test");
julia> println("Section '"* s.metadata["name"] * "' has entries with salinities starting:")
Section 'Test' has entries with salinities starting:
julia> for i in 1:length(s.data)
           println(first(s.data[i].data.salinity,3))
       end
[29.921, 29.9205, 29.9206]
[30.921, 30.9205, 30.9206]
source
OceanAnalysis.beam_to_xyzFunction
beam_to_xyz(adp::Adp; debug::Int64=0)

Change velocity in an RDI Adp object from beam to xyz coordinates

This is done by using the transformation_matrix that is stored within adp. See read_adp_rdi for how to read an RDI Adp object, and xyz_to_enu for how to convert it from xyz to enu coordinates.

Examples

using OceanAnalysis, Plots
file = joinpath(dirname(dirname(pathof(OceanAnalysis))),
    "data", "adp_rdi.000")
adp = read_adp_rdi(file);
adp_xyz = beam_to_xyz(adp);
v = adp.data["velocity"];
V = adp_xyz.data["velocity"];
CM = cgrad(:RdBu, rev=true);
p1 = heatmap(transpose(v[:, :, 1]), c=CM, title="beam 1", titlefontsize=9);
p2 = heatmap(transpose(v[:, :, 2]), c=CM, title="beam 2", titlefontsize=9);
p3 = heatmap(transpose(v[:, :, 3]), c=CM, title="beam 3", titlefontsize=9);
p4 = heatmap(transpose(v[:, :, 4]), c=CM, title="beam 4", titlefontsize=9);
pu = heatmap(transpose(V[:, :, 1]), c=CM, title="u", titlefontsize=9);
pv = heatmap(transpose(V[:, :, 2]), c=CM, title="v", titlefontsize=9);
pw = heatmap(transpose(V[:, :, 3]), c=CM, title="w", titlefontsize=9);
pe = heatmap(transpose(V[:, :, 4]), c=CM, title="err", titlefontsize=9);
plot(p1, p2, p3, p4, layout=(4, 1), size=(1000, 700))
plot(pu, pv, pw, pe, layout=(4, 1), size=(1000, 700))
source
OceanAnalysis.coastlineFunction
coastline(symbol::Symbol=:world_fine)

Access a built-in Coastline dataset. The only valid choices for name are :global_coarse and :global_fine. These are handled by reading the built-in datasets data/coastline_coarse.csv.gz and datasets data/coastline_fine.csv.gz, respectively.

Examples

# Nova Scotia
using OceanAnalysis, Plots
cl = coastline(:global_fine);
plot_coastline(cl, xlims=(-68, -58), ylims=(43, 48))
source
coastline(filename::String, header::Integer=0)

Create a Coastline object from a CSV file.

The file must have columns named longitude and latitude, with NaN values indicating breaks separating islands, etc.

The work is done by passing filename and header to CSV.read(). The file must have 1 or more header lines, the last of which must contain column names longitude and latitude. NaN values will be taken to indicate breaks between segments that trace continents, nations, etc.

Examples

# World view
using OceanAnalysis, Plots
dir = dirname(dirname(pathof(OceanAnalysis)))
file = joinpath(dir, "data", "coastline_coarse.csv.gz")
cl = coastline(file, 1)
plot_coastline(cl)
source
coastline(longitude::Union{AbstractVector,AbstractRange},
    latitude::Union{AbstractVector,AbstractRange})

Create a Coastline object from longitude and latitude values.

Use NaN values for both longitude and latitude to indicate breaks in the coastline from continent to continent, nation to nation, etc.

Examples

# Nova Scotia
using OceanAnalysis, CSV, Plots, DataFrames
dir = dirname(dirname(pathof(OceanAnalysis)));
file = joinpath(dir, "data", "coastline_fine.csv.gz");
data = CSV.read(file, DataFrame, header=1);
cl = coastline(data.longitude, data.latitude);
plot_coastline(cl, xlims=(-68, -58), ylims=(43, 48))
source
OceanAnalysis.coriolisFunction
coriolis(latitude::Real)

Compute the value of the Coriolis parameter, given a latitude in degrees.

Examples

julia> coriolis(45.0)
0.00010312607931384281
source
OceanAnalysis.coordinate_from_stringFunction
degree = coordinate_from_string(s::String)

Try to extract a longitude or latitude from a string. If there are two (space-separated) tokens in the string, the first is taken as the decimal degrees, and the second as decimal minutes. The goal is to parse hand-entered strings, which might contain letters like "W" and "S" (or the same in lower case) to indicate the hemisphere. Humans are quite good at writing confusing strings, so this function is not always helpful.

Examples

julia> using OceanAnalysis

julia> coordinate_from_string("1.5")
1.5

julia> coordinate_from_string("1 30")
1.5

julia> coordinate_from_string("1S 30")
-1.5

julia> coordinate_from_string("27* 14.072 N")
27.234533333333335

julia> coordinate_from_string("111* 31.440 W")
-111.524
source
OceanAnalysis.CTFunction
CT(SA, temperature, pressure)
CT(ctd)

Compute Conservative Temperature (CT), using gsw_ct_from_t() in the GibbsSeaWater package.

The first form takes single values and returns a single value. The second form extracts values from a Ctd object and then calls the first form as CT.() so that it returns a vector of CT values.

Examples

julia> using OceanAnalysis

julia> CT(35.0, 10.0, 100.0)
9.981322531922249
source
OceanAnalysis.despikeFunction
despike(x; k::Int64=7, n::Int64=4, action::Symbol=:replace)

Despike a timeseries, or reveal spikes.

The procedure starts by computing a running median version (with window length k) of the provided time-series. Then the standard deviation of the difference between the two timeseries is computed. Points are considered to be spikes if they have difference exceeding n times the overall standard deviation. Any such points are then replaced with the running-median values, and the resultant possibly-altered timeseries is returned. The default values of k and n may provide a good starting point, but users are advised to explore other values, in the context of the character of data being analyzed.

Parameters

  • x a vector of values representing a timeseries.

Keywords

  • k width of running-median filter. This is passed to running_median, which computes the running mean. Since the computation uses centred windows, k will be incremented by 1 if it is an even number.
  • n spike criterion. Increasing n will usually decrease the number of points considered spikes, unless they are very anomalous.
  • action a Symbol indicating what to do with spikes. If this is :replace (which is the default) then spikes are replaced with the corresponding running-median values. If it is :NaN then they are replaced with NaN values. And, if it is :reveal then a BitVector is returned, with 0 at the indices of non-spike values and 1 at the indices of spike values.

Return value

Either a Float64 vector (if action is :replace or :NaN) or a logical BitVector (if action is :reveal).

Examples

using OceanAnalysis, Plots
i = 1:40;
x0 = sin.(2 * pi * i / 40);
x = x0 .+ (rand(length(x0)) .- 0.5) / 10.0;
x[10] = x[10] + 1;
xd = despike(x);
# Plot 'base' (before noise), 'signal' (base + noise) and 'despiked'
plot(i, x0, label="base", ms=3)
scatter!(i, x, label="signal")
scatter!(i, xd, label="despiked", ms=2)
# Print overview of spikes (show it and nearest neighbours)
spike_indices = findall(despike(x, action=:flag));
for index in spike_indices
    println(x[index.+range(-1, 1)])
end
source
OceanAnalysis.geod_distanceFunction
geod_distance(lon1, lat1, lon2, lat2; a=6378137.00, f=1.0 / 298.257223563)

Return the geodetic distance between two points on earth, in kilometres.

This uses the Vincenty [1975] formulation. According to tests described in that paper, the results can be expected to be accurate to 0.01 mm.

Note that geod_distance() mimics the geodDist() function in the oce R package (Ref. 2), and is based on the R and C++ code used in the latter. The code trail is 5 decades long, starting with a CDC-6600 Fortran version written in 1974 by L. Pfeifer, which was modified for IBM 360 by J. G. Gergen. In 2003, D. Gillis wrote an R version that was modified for inclusion in the oce package by Dan Kelley.

Arguments

  • lon1 and lat1, in degrees, give the location of a given point on the earth.
  • lon2 and lat2, in degrees, give the location of a given point on the earth.
  • a and f are the semi-major axis and flattening parameter of the reference ellipsoid. (The default values are highly recommended.)

References

  1. Vincenty,T. 1975. Direct and inverse solutions of geodesics on the ellipsoid with application of nested equations. Survey Review 23(176):88-94. https://www.ngs.noaa.gov/PUBS_LIB/inverse.pdf

  2. Kelley, Dan E., Clark Richards, and Chantelle Layton. “Oce: An R Package for Oceanographic Analysis.” Journal of Open Source Software 7, no. 71 (2022). https://doi.org/10.21105/joss.03594

source
OceanAnalysis.get_amsrFunction
get_amsr(date::String)
source
get_amsr(date::Date=Dates.today() - Day(4); type::String="3day",
    destdir::String=".", server::String="https://data.remss.com/amsr2/ocean/L3/v08.2",
    debug::Integer=0)

Download a Advanced Microwave Scanning Radiometer data file.

This works by constructing a filename to be downloaded. If that file does not exist in destdir, then it is downloaded from the server, and get_amsr returns the full path to that existing file. Otherwise, the file is downloaded, and the return value is the path to the resultant local file.

For example, if this function is called on 2025-09-27 with no arguments specified, an attempt will be made to download a file named "./RSS_AMSR2_ocean_L3_3day_2025-09-23_v08.2.nc" from the server `https://data.remss.com/amsr2/ocean/L3/v08.2/3day/2023/".

See read_amsr for how to deal with the files downloaded by get_amsr.

Arguments

  • date a Date object specifying the requested measurement time. This defaults to 4 days prior to the current date.

Keywords

  • destdir: Path to the destination directory. The author usually sets destdir="~/data/amsr", so that the file will be in a central location for use by other analysis procedures.

  • server: The base of the server location. The default value ought to be used unless the data provider changes their web scheme, although the likelihood of the query working in such a case is slim, since changes tend to be sweeping.

  • type: The type of data requested. At the moment, the only choice is "3day" (the default), for a composite covering 3 days of observation, which removes most viewing-path and cloud blanks. If there is sufficient need, other types may be added, from the list: "daily" for a daily reading, "weekly" for a composite covering a week, and "monthly" for a composite covering a month. In the "daily" case, the data arrays are 3D, with the third dimension representing ascending and descending traces, but in all the other cases, the arrays are 2D.

  • debug: An indication of whether to print information during processing. The default value of 0 means to work quietly, and any larger integer indicates to print some information.

Return

get_amsr returns a string that is the full pathname of the downloaded file, which may be supplied as the first argument to a call to read_amsr.

source
OceanAnalysis.get_argo_indexFunction
get_argo_index(destdir::String="."; age::Real=1.0,
    server::String="https://data-argo.ifremer.fr", debug::Int64=0)

Download an Argo index file, unless an existing local copy is newly downloaded.

The file is obtained from server and stored in destdir, but only if an existing version of the file in destdir is under age days old. The default is to cache index files for one day.

The default destdir is the local directory, but it is common to set this value to some central location so the file can be used by Julia sessions running in multiple directories. (The author sets destdir="~/data/argo", for example.)

Return value

read_argo_index returns the full name of the downloaded (or cached) file.

Use read_argo_index to interpret the downloaded file.

source
OceanAnalysis.get_argoFunction
get_argo(filename::String="", destdir::String="."; age::Real=1.0; server::String="https://data-argo.ifremer.fr", debug=0)

Download an Argo profile file, if an existing copy is less than age days old.

Arguments

  • file: path to the file on a server, of the form agency/#/profiles/-#_##.nc, where # is the ID number of the float, - tells the status of the file (R for realtime files, D for delayed-mode files, etc.) and ## is an identifier for the cast (usually but not always an integer value). This system matches the file information stored on argo servers, as downloaded by get_argo_index and read by read_argo_index.

Keywords

  • destdir: name of the directory into which to save the file.

  • age: file-caching time in days. If the requested file does not exist locally then age is ignored and the file is downloaded. It will also be downloaded if there is an existing file but it was last downloaded more than age days ago.

  • server: base URL for the server (it is very unlikely that users will set this

  • debug: integer indicating whether to print information during processing. The default

value of 0 means to work silently.

Returns

  • get_argo returns the full path name of the local file after downloading, or as cached recently.

Example

# Get most last-named file in the Argo index and plot a temperature profile
using OceanAnalysis
index_file = get_argo_index()
index = read_argo_index(index_file)
argo_file = get_argo(index.file[end])
argo = read_argo(argo_file)
plot_profile(argo, which="CT")
source
OceanAnalysis.get_elementFunction
get_element(x::OA, element::String; debug)

Get an element from an object.

If x.metadata holds an item of the given name, return that item. If not, and if x.data holds such an item, return that item.

If x is a Ctd object, some computed things may be returned. These are N2, SA, CT, sigma0 and spiciness0.

If x is a Section object, then get_element can return any item from the metadata of the constituent Ctd objects that are stored in x.data.

source
OceanAnalysis.get_fileFunction
get_file(url::String=""; destdir::String=".", age::Real=1.0, debug::Int64=0)

Download/cache a remote file

The url is a remote source for a file to be downloaded. If no file of that name exists in the destination directory, destdir, then the file is downloaded and its name is returned. However, if such a file already exists, and if its age is under age days, then the file is assumed to be up-to-date and is not downloaded. Some processing steps are printed if debug>0.

source
OceanAnalysis.get_nc_valueFunction
get_nc_value(nc, name; require_valid=false)

Transform an item from a NetCDF file into a more useable object.

Any missing values are converted to NaN. This function is designed only for numerical values, not strings or other types.

Arguments

  • nc a value returned by NCDataset().

  • name the name of an object contained in nc.

Keywords

  • require_value boolean value indicating whether to report an error if the desired element consists entirely of bad values.
source
OceanAnalysis.get_sectionFunction
get_section(url::String; destdir=".", debug::Bool=false)

Download a zipfile from url and expand its contents into files within a destination directory inferred from the URL or as defined by destdir, if the latter is not missing.

The work starts by downloading a zipfile to the local directory, if it is not already present. Then a directory name is constructed based on url and the value of destdir. If no such directory exists, it is created. Then the zipfile is expanded, storing the contents in this new directory.

The return value is the name of the new directory.

Setting debug=1 will cause the printing of some of the processing steps.

Examples

using OceanAnalysis
# Saves files into a local directory called 'ar07_74JC20140606'.
url = "https://cchdo.ucsd.edu/data/11852/ar07_74JC20140606_ct1.zip"
sdir = get_section(url)
println("Downloaded ", length(readdir(sdir)), " files to '", sdir, "'")
source
OceanAnalysis.get_topographyFunction
get_topography(name::Symbol=:global_coarse; debug::Int64=0)

Access a built-in topography dataset.

The only valid choice for name is :global_coarse, which returns a world view at resolution of 30 minutes, i.e. with grid cells near the equator spanning approximately 30 nautical miles.

See also read_topography.

source
get_topography(west::Real, east::Real,
    south::Real, north::Real; resolution::Real=4.0, destdir::String = ".",
    server::String = "https://gis.ngdc.noaa.gov", debug::Int64 = 0)

Download and cache a topography file.

Topographic data are downloaded from a data server that holds the ETOPO1 dataset (see Amante and Eakins, 2009, for an introduction to the data and see Pante and Simon-Bouhet, 2013, for code that queries a server in a manner similar to that used here), and saved as a netCDF file that has a name that reveals the data request, if a file of that name is not already present on the local file system. The return value is the name of the data file, and its typical use is as the filename for a call to read_topography. Subsequent calls to get_topography with identical parameters will return the name of an already-downloaded file, without downloading a new copy.

The specified longitude and latitude limits are rounded to 2 digits after the decimal place (corresponding to an equatorial footprint of approximately 1 km), and these are used in the server request.

The region of interest is defined by a rectangle bounded by the values of west, east, south and north. The resolution is set by the value of resolution, which is minutes. The default, 4.0 minutes, corresponds to 4 nautical miles (approx. 7.4km) in the north-south direction, and less in the east-west direction. The default values of these 5 arguments yield a view of the Bay of Fundy region.

The value of server ought not to be modified by users, except perhaps to experiment if the NOAA server changes. (It is unlikely that merely changing this value will help much, though, since changes tend not to be small on these servers.)

References

  1. Amante, C. and B.W. Eakins, 2009. ETOPO1 1 Arc-Minute Global Relief Model:

Procedures, Data Sources and Analysis. NOAA Technical Memorandum NESDIS NGDC-24. National Geophysical Data Center, NOAA. doi:10.7289/V5C8276M

  1. Pante, Eric, and Benoit Simon-Bouhet. "Marmap: A Package for Importing,

Plotting and Analyzing Bathymetric and Topographic Data in R." PLoS ONE 8, no. 9 (2013): e73051. doi:10.1371/journal.pone.0073051. (The package referenced was updated on 2025-Aug-2; for the query generation, see the fetch function of that package's source code in R/getNOAA.bathy.

  1. API https://gis.ngdc.noaa.gov/arcgis/help/en/rest/services-reference/enterprise/export-image/
source
OceanAnalysis.gravityFunction
gravity(latitude::Real=45.0, z::Real=0.0)

Gravitational acceleration for an ellipsoidal model of the earth.

Computes g (in m/s^2) based on latitude (in degrees North) and elevation z (in m above the geoid), using the same formula as that in the first form of equation D.3 in reference 1.

References

  1. Intergovernmental Oceanographic Commission, SCOR, and IAPSO. The International Thermodynamic Equation of Seawater–2010: Calculation and Use of Thermodynamic Properties. No. 56. UNESCO, 2016.

Examples


julia> gravity()
9.806199977310339
source
OceanAnalysis.grid_ctdFunction
grid_ctd(ctd::Ctd;
    pressure_grid::Vector{Float64}=missing, pressure_step::Float64=2.0;
    method::Symbol=:interpolate, debug::Int64=0)

Grid a CTD to standardized pressure levels.

The levels are as given by pressure_grid or, if that is not provided, as the pressure range within ctd, incrementing by pressure_step.

Arguments

  • ctd a Ctd to be gridded. It must contain a pressure column.

Keywords

  • method: a symbol indicating the gridding method. At the moment, only one choice is accepted, namely :interpolate, which means to use linear interpolation of each field to the specified pressure grid.

  • debug: an optional value that, if it exceeds 0, indicates that debugging output should be printed during processing.

Example

using OceanAnalysis, Plots
f = joinpath(dirname(dirname(pathof(OceanAnalysis))), "data", "ctd.cnv");
ctd = read_ctd_cnv(f);
pressure_grid = 0.0:1.0:100.0;
ctd2 = grid_ctd(ctd, pressure_grid=pressure_grid);
plot(ctd["salinity"], -ctd["pressure"], seriestype=:path, legend=false, framestyle=:box)
plot!(ctd2["salinity"], -ctd2["pressure"], seriestype=:path, color=:red, legend=false, framestyle=:box)
source
OceanAnalysis.grid_sectionFunction
grid_section(section::Section, pressure_step::Float64=2.0; debug::Int64=0)

Grid a section, altering all the CTD objects to employ a uniform pressure grid.

The pressure grid is set up to range from 0 to the maximum pressure across all the CTDs in section, incrementing by pressure_step. Once this is set up, the other fields are interpolated to this grid using linear interpolation, with missing values inserted for gridded pressures that exceed the maximum value in any given CTD.

Arguments

  • section a Section object.

  • pressure_step desired pressure increment, in dbar.

Keywords

  • debug: an optional value that, if it exceeds 0, indicates that debugging output should be printed during processing.

Examples

using OceanAnalysis
url = "https://cchdo.ucsd.edu/data/11852/ar07_74JC20140606_ct1.zip"
dir = get_section(url)
section = read_section(dir);
section_gridded = grid_section(section)
source
OceanAnalysis.interpolate_barnesFunction
interpolate_barnes(
    x::Union{AbstractVector,AbstractRange},
    y::Union{AbstractVector,AbstractRange},
    z::Union{AbstractVector,AbstractRange},
    w::Union{AbstractVector,AbstractRange,Nothing}=nothing;
    xg::Union{AbstractVector,AbstractRange,Nothing}=nothing,
    yg::Union{AbstractVector,AbstractRange,Nothing}=nothing,
    xr::Union{Float64,Nothing}=nothing,
    yr::Union{Float64,Nothing}=nothing,
    gamma::Float64=0.5, iterations::Int64=2, debug=0)

Interpolate a two-dimensional field to a grid, using the Barnes method as described by Barnes (1994a, 1994b, 1994c). The methodology follows that of the interp_barnes function in the R oce package.

Arguments

  • x Vector of Float64 values for the x coordinate of data points.
  • y Vector of Float64 values for the z coordinate of data points.
  • z Vector of Float64 values for the z values at the (x,y) data points.
  • w Vector of Float64 values for the weights for the data points. If not provided, this defaults to a vector of unit values.

Keywords

  • xg Float64 vector giving grid coordinates in the x direction. If not provided, this is computed with pretty(x,50).
  • yg As xg but for the y direction.
  • xr Float64 telling the influence scale in the x direction. If not provided, this defaults to the range of x values, divided by the square root of the number of x values.
  • yr As xr but for the y direction.
  • gamma Float64 telling the value of gamma to use (0.5 by default).
  • iterations integer telling how many iterations to perform (2 by default).
  • debug an integer indicating whether to print information during processing. The default value of 0 means to work quietly, and any larger integer indicates to print some information.

Value

A Dict with elements named xg, yg and zg that hold the grid coordinates and values, along with wg (the final weights) and zd (the values interpolated at the data coordinates).

Examples

using OceanAnalysis, CSV, DataFrames, Statistics, Plots
file = joinpath(dirname(dirname(pathof(OceanAnalysis))), "data", "wind.csv")
data = CSV.read(file, DataFrame);
w = repeat([1.0], nrow(data));
xg = range(0.0, 11.0, step=0.2);
yg = range(0.0, 9.0, step=0.2);
res = interpolate_barnes(data.x, data.y, data.z)
scatter(data.x, data.y, framestyle=:box, label=false, ms=2, tickdirection=:out,
    xlab="x", ylab="y", xlim=(0, 11), ylim=(0, 9))
annotate!(data.x, data.y .+ 0.2, text.(data.z, 8, :blue))
contour!(res["xg"], res["yg"], res["zg"],
    levels=10:5:30, cbar=false, clabels=true, c=:black)

References

Barnes, Stanley L. “Applications of the Barnes Objective Analysis Scheme. Part I: Effects of Undersampling, Wave Position, and Station Randomness.” Journal of Atmospheric and Oceanic Technology. Journal of Atmospheric and Oceanic Technology 11, no. 6 (1994a): 1433–48.

Barnes, Stanley L. “Applications of the Barnes Objective Analysis Scheme. Part II: Improving Derivative Estimates.” Journal of Atmospheric and Oceanic Technology. Journal of Atmospheric and Oceanic Technology 11, no. 6 (1994b): 1449-1458.

Barnes, Stanley L. “Applications of the Barnes Objective Analysis Scheme. Part III: Tuning for Minimum Error.” Journal of Atmospheric and Oceanic Technology. Journal of Atmospheric and Oceanic Technology 11, no. 6 (1994c): 1459–79.

source
OceanAnalysis.label_from_varnameFunction
label_from_varname(String::varname, abbreviate::Symbol=:long)

Return a string suitable for use as an axis label.

Parameters

  • varname String holding the variable name.
  • abbreviate Symbol indicating the available size, with valid choices being :short, :medium and :large.

Examples

using OceanAnalysis
label_from_varname("CT", :short)
label_from_varname("CT", :medium)
label_from_varname("CT", :large)
label_from_varname("CT") # defaults to :large
source
OceanAnalysis.N2Function
N2(ctd::Ctd, s::Float64=0.15; debug::Int64=0)

Compute the square of the buoyancy frequency, N², for a Ctd object. The value is inferred from a smoothing cubic spline that models the dependence of sigma0 on pressure.

In the present version, the spline is fitted with the Dierckx::Spline1D() function (Reference 1), which is provided with equal weights, w, for all points, with k=3 to set the polynomial order to cubic, and with the user-specified bc to control behaviour near top and bottom, along with s to control smoothness. After the derivative is computed, any negative N² values are set to zero.

Parameters

  • ctd a [Ctd] object

Keywords

  • s either a Float64 value or a symbol. In the first case, it is the value of s supplied to Dierckx::Spline1D(), which is used to smooth the density curve as a function of pressure. According to the documentation for the Fortran code behind this function (see https://www.netlib.org/dierckx/curfit.f), a reasonable starting point for exploring the dependence of the spline curve on s is in the range from delta^2*(N-sqrt(2N)) to delta^2*(N+sqrt(2N)), where N is the number of data points and delta is a measure of the density "wiggles" across which to smooth the spline. This is the value used if s=:auto (the default) is specified; using s=:smooth and s=:rough multiplies this value by sqrt(2) and 1/sqrt(2), respectively. Note that specifying debug=1 will make the function print out the numeric value of s used, if one of these three symbols has been supplied.

  • bc a string, passed to Dierckx::Spline1D(), that indicates what to do near boundaries. The default is "nearest".

  • debug an integer indicating whether to print information during processing. The default value of 0 means to work quietly, and any larger integer indicates to print some information.

References

  1. https://github.com/JuliaMath/Dierckx.jl

Examples

# Demonstrate N2()
using OceanAnalysis, Plots
pkgdir = dirname(dirname(pathof(OceanAnalysis)))
filename = joinpath(pkgdir, "data", "ctd.cnv")
ctd = read_ctd_cnv(filename);
histogram(N2(ctd), label="N²")
plot(p1, p2)
source
OceanAnalysis.OceanAnalysisModule

The OceanAnalysis module is intended to help with the analysis of oceanographic data. It is in a preliminary form, providing help with only four file types: (1) Argo NetCDF files, (2) CTD files in the Seabird CNV format, (3) AMSR satellite NetCDF files and (4) coastline files. In each case, the capabilities are quite limited, reflecting the early stage of the package. Users who need more powerful tools for reading and analysing oceanographic data, ought consider using the oce package in the R language, which has more capabilities and has been in wide use for over a decade.

The functions that read data return objects that are structs holding two items: data and metadata.

  • The data item holds the actual data. The form of the data depends on the

class of the object. For example, in a Ctd object, data is a DataFrame, whereas in a Section object, data is a Vector of Ctd objects.

  • The metadata item holds information about the data. For example, with a

Ctd object, metadata holds the location and time of sampling, along with other information, depending on the source of the data.

As a convenience, get_element may be used to extract information from either the metadata or data parts of an OceanAnalysis object.

FIXME: discuss the fact that get_element can return derived quantities.

FIXME: discuss [ here. I think maybe we ought to cause [ to call get_element regardless of the object class. But, in that case, is there any need for get_element at all?

source
OceanAnalysis.plot_adpFunction
plot_adp(adp::Adp; which=:velocities, debug::Int64=0, kwargs...)

Plot the data stored in an Adp object.

This function provides some basic plots of the contents of an acoustic-Doppler profiler (Adp) object.

Arguments

  • adp an Adp, as created with read_adp_rdi.

  • which a Symbol indicating what to plot. If which is :velocity1 then a [heatmap] plot is made of the first component of velocity. It will be labelled as "beam 1", "ũ" or "u", according to whether adp["coordinate_system"] is :beam, :xyz or :enu. Similar results are obtained for :velocity2 etc., where the fourth element is called "ẽ" or "e", designating an error estimate. If which is velocities, then the result is a multi-panel plot, with one panel per velocity component. If which is :heading then a time-series plot is made of heading, with analogous results for :pitch and :roll. If which is :angles then a three-panel plot is made, showing these three angles. If which is :uv and adp["coordinate_system"] is :enu, then mid-distance east and north components of velocity are computed and then plotted in a scatterplot.

Keywords

  • debug: an optional integer value that, if it exceeds 0, indicates that debugging output should be printed during processing.

  • kwargs: optional items, passed to heatmap for velocity fields, or to scatter for time-series and other x-y plots.

Examples

using OceanAnalysis, Plots
file = joinpath(dirname(dirname(pathof(OceanAnalysis))),
    "data", "adp_rdi.000")
adp_beam = read_adp_rdi(file);
plot_adp(adp_beam)
adp_xyz = beam_to_xyz(adp_beam);
plot_adp(adp_xyz)
source
OceanAnalysis.plot_amsrFunction
plot_amsr(amsr::Amsr;
    xlims=[0.0, 360.0], ylims=[-90.0, 90.0], tickdirection=:out,
    color=:turbo, levels=[], clim=:auto, size=(800, 550), dpi=300,
    debug::Int64=0)

Plot a heatmap of an AMSR field. By default, SST is shown using the default Julia colorscheme, and the view is of the whole earth. See the example for how to use another colorscheme, and how to narrow the geographical view. Note that the graph scales longitude and latitude so that a hypothetical circular island at the midpoint of the view would be drawn as a circle.

Arguments

Keywords

  • xlims: The range of longitude to be shown. This is based on the 0-to-360 notation, since that is how AMSR data are stored.

  • ylims: The range of latitude to be shown.

  • tickdirection: The direction of axis tick marks. The default is for them to point outward, opposite to the Julia default.

  • color: The colour scheme for the heatmap. The default, :turbo, is a rainbow-like scheme. Other popular choices include :viridis for a green-hued scheme, and :auto for the default yellow-hued Julia scheme.

  • levels: A vector holding the desired contour levels, a single integer giving the desired number of contours, or :auto for automatic selection. The default is useful for world views of sea-surface temperature.

  • clim: A tuple specifying the range of values to be represented by the color scheme. If not provided, this defaults to the range of the data in the chosen view.

  • size: A numeric tuple holding the size of the plot.

  • dpi: A number representing the resolution of the plot, in dots per inch.

  • debug: An integer controlling whether to print information during processing. The default is to work silently; use any positive value to get some printing.

Examples

using OceanAnalysis
file = "~/data/amsr/RSS_AMSR2_ocean_L3_3day_2025-09-07_v08.2.nc"
amsr = read_amsr(file, "SST");
plot_amsr(amsr, xlims=(300,360), ylims=(40,60))
source
OceanAnalysis.plot_coastlineFunction
plot_coastline(coastline::Coastline; debug::Int64=0, kwargs...)

Plot a coastline with longitude and latitude axes (i.e. without a map projection).

The aspect_ratio of the plot is set as the reciprocal of the mean of the ylims values, to preserve shapes near that spot.

Arguments

Keywords

  • debug an integer indicating whether to print information during processing. The default value of 0 means to work quietly, and any larger integer indicates to print some information.

  • kwargs... other arguments, passed to plot, e.g. xlim and ylim to control the plot view, color for the land colour, etc.

Examples

using OceanAnalysis, Plots
#plot_stations([-59.9], [43.9]; xlim=(-70, -50), ylim=(40, 50), color=:green, markercolor=:red)
# Default world view
world = plot_coastline(coastline())
# Nova Scotia view, with light-gray for land
ns = plot_coastline(coastline(); color=:gray80, xlim=(-70.0, -55.0), ylim=(43.0, 48.0))
plot(world, ns)
source
OceanAnalysis.plot_coastline!Function
plot_coastline!(coastline::Coastline;
    seriestype=:shape, color=:bisque3, linewidth::Real=0.5)

Add a coastline to an existing plot.

This shares several arguments with plot_coastline, but not those that could alter the geometry. Note that the plot limits are inherited from the existing plot, so xlim and ylim should not be supplied in the kwargs... grouping.

Arguments

Keywords

  • debug an integer indicating whether to print information during processing. The default value of 0 means to work quietly, and any larger integer indicates to print some information.

  • kwargs... other arguments, passed to plot, e.g. xlim and ylim to control the plot view, color for the land colour, etc.

source
OceanAnalysis.plot_echosounderFunction
plot_echosounder(e::Echosounder; debug::Int64=0, kwargs...)

Plot the data stored in an Echosounder object.

This function provides some basic plots of the contents of an acoustic-Doppler profiler (Adp) object.

Arguments

  • e an Echosounder object, as created with read_echosounder.

  • which a Symbol indicating what to plot. If which is :log_amplitude then a [heatmap] plot is made of the base-10 logarithm of the signal amplitude. Eventually, other possible values of which may be handled.

Keywords

  • debug: an optional integer value that, if it exceeds 0, indicates that debugging output should be printed during processing.

  • kwargs: optional items, passed to heatmap.

Examples

using OceanAnalysis
f = "/Users/kelley/Dropbox/data/archive/sleiwex/2008/fielddata/2008-07-01/Merlu/Biosonics/20080701_163942.dt4";
if isfile(f)
    e = read_echosounder(f);
    plot_echosounder(e)
end
source
OceanAnalysis.plot_profileFunction
plot_profile(ctd::Ctd; which::String="CT", vertical::Symbol=:pressure,
    abbreviate::Symbol=:long, debug::Int64=0, kwargs...)

Plot an oceanographic profile for data contained in ctd, showing how the variable named by which depends on either pressure or density. The variable is drawn on the x axis and pressure on the y axis. Following oceanographic convention, the y axis is set up so that waters nearer the air-sea interface are nearer the top of the plot.

Arguments

  • ctd a Ctd object to be plotted.

Keywords

  • which an indication of what to plot on the x axis. The default value, "CT", indicates to plot Conservative Temperature. Anything stored in the object's data can be plotted, along with some things that can be calculated from these values. Common choices include: "N2" for N², the square of the buoyancy frequency; "SA" for the Gibbs Seawater formulation of Absolute Salinity; "salinity" for Practical Salinity; "sigma0" for the Gibbs Seawater formulation of density anomaly referenced to the surface; "spiciness0" for the Gibbs Seawater seawater spiciness referenced to the surface; and "temperature" for in-situ temperature.

  • vertical a Symbol specifying what to plot on the y axis. The default is :pressure, but :density is also permitted.

  • abbreviate a Symbol indicating a category for axis length, used in determining how to label the axes. The valid choices are :short, :medium, and :long.

  • debug indicator of debugging level. If this exceeds 0, some information is printed during processing.

  • kwargs... is passed to plot(), to permit further customization; see https://docs.juliaplots.org/stable/ for more information on possibilities.

Examples

using OceanAnalysis, Plots

# Example 1: show overview of an Argo profile
pkgdir = dirname(dirname(pathof(OceanAnalysis)))
f = joinpath(pkgdir, "data", "D4902911_095.nc")
ctd = read_argo(f) |> as_ctd;
# Plot profiles of Conservative Temperature, Absolute Salinity, and potential
# density anomaly with respect to surface pressure.
p1 = plot_profile(ctd; which="CT")
p2 = plot_profile(ctd; which="SA")
p3 = plot_profile(ctd; which="sigma0")
plot(p1, p2, p3, layout=(1, 3), size=(800, 400))

# Example 2: add a new variable to the profile, then plot it
using GibbsSeaWater
ctd.data.conductivity = gsw_c_from_sp.(ctd["salinity"], ctd["temperature"], ctd["pressure"]);
plot_profile(ctd, which="conductivity", xlab="Conductivity [mS/cm]")
source
OceanAnalysis.plot_sectionFunction
plot_section(section::Section, which="salinity";
    type=:contourf, xvar=:latitude, yvar=:pressure, debug::Int64=0, kwargs...)

Arguments

  • section a Section, as created with as_section or read_section.

  • which a String indicating the name of the hydrographic variable to be plotted. This must be present in each of the Ctd objects stored within the section.data. Another requirement is that the section has been gridded, using grid_section. The plotting is done with contour, contourf or heatmap as directed by the type argument. In each case, kwargs... is passed to the function to permit customization. If show_stations is true, then vline is used to draw vertical lines indicating station locations. Note that case 2, section_is_gridded is called first to ensure that the section has been gridded with grid_section, with an error being reported if not.

Keywords

  • type a Symbol indicating the type of plot. This may be :contour for simple contours, :contourf (the default) for filled contours, or :heatmap for an image.

  • xvar either a Symbol (which must be one of :distance, :latitude or :longitude) or a Tuple with two elements the first being a label for the x axis and the second being a vector of values for x that correspond to the stations in section.data. See the Examples for a case in which sampling time is used for the Tuple case.

  • yvar a Symbol, the permitted values of which are :depth and :pressure.

  • show_stations a Bool value indicating whether to draw vertical gray dotted lines to indicate station locations on cross-section diagrams.

  • debug: an optional integer value that, if it exceeds 0, indicates that debugging output should be printed during processing.

  • kwargs: optional items, passed down to lower-level plotting functions. For example, size controls the size of the plot, xlim and ylim control the viewing window, and color controls the colour.

Examples

using OceanAnalysis, Plots
url = "https://cchdo.ucsd.edu/data/41926/90CT40_1_ct1.zip"; # exchange format
dir = get_section(url);
s = read_section(dir);
s.data = s.data[s["longitude"].< (-68.0)];
p1 = plot_stations(s, xlim=(-80, -65), ylim=(35, 43));
scale_bar(500);
# Note that we must grid to get the cross-section diagrams
sg = grid_section(s);
p2 = plot_section(sg, "salinity", ylim=(0, 2000));
p3 = plot_section(sg, "salinity", ylim=(0, 2000), xvar=("Time", sg["time"]));
l = @layout [a; b c];
plot(p1, p2, p3, layout=l, dpi=300, size=(800, 500))
source
OceanAnalysis.plot_stationsFunction
plot_stations(longitude::Vector{Float64}, latitude::Vector{Float64};
    draw_coastline::Bool=true, debug::Int64=0, kwargs...)

Plot station locations (or similar data) on a map.

The locations are plotted with scatter, with aspect ratio computed using the midpoint between the extrema of latitude. If show_coastline is true, then plot_coastline is called to draw the land.

Arguments

  • longitude a vector of longitudes, in the -180 to 180 degree frame.

  • latitude a vector of latitudes, in the -90 to 90 degree frame.

Keywords

  • draw_coastline: a Bool value indicating whether to plot the coastline, using plot_coastline.

  • debug: an optional value that, if it exceeds 0, indicates that debugging output should be printed during processing.

  • kwargs: optional items, passed down to lower-level plotting functions. For example, size controls the size of the plot, xlim and ylim control the viewing window, color controls the land colour, and markercolor controls the station-location colour.

Examples

using OceanAnalysis, Plots
# Twenty fake stations between Halifax and Sable Island
lon = collect(range(-59.91, -63.53, length=20))
lat = collect(range(43.93, 44.59, length=20))
plot_stations(lon, lat, xlim=(-65.0, -59.0), ylim=(43.0, 46.0))
source
plot_stations(section::Section; draw_coastline::Bool=true, debug::Int64=0, kwargs...)

Plot section station locations on a map.

The locations are plotted with scatter, with aspect ratio computed using the midpoint between the extrema of latitude. If show_coastline is true, then plot_coastline is called to draw the land.

Arguments

Keywords

  • draw_coastline: a Bool value indicating whether to plot the coastline, using plot_coastline.

  • debug: an optional value that, if it exceeds 0, indicates that debugging output should be printed during processing.

  • kwargs: optional items, passed down to lower-level plotting functions. For example, size controls the size of the plot, xlim and ylim control the viewing window, color controls the land colour, and markercolor controls the station-location colour.

Examples

using OceanAnalysis, Plots
url = "https://cchdo.ucsd.edu/data/41926/90CT40_1_ct1.zip"; # exchange format
dir = get_section(url);
s = read_section(dir);
plot_stations(s, xlim=(-80,0), ylim=(20,50))
source
OceanAnalysis.plot_TSFunction
plot_TS(ctd::Ctd; sigma0_levels=[], spiciness0_levels=0,
    draw_freezing=true, abbreviate=false, debug::Int64=0, kwargs...)

Plot an oceanographic TS diagram, with the Gibbs Seawater equation of state.

Whether contours of density and spiciness are drawn depends on values of the sigma0_levels and spiciness0_level arguments. There are 4 categories. (1) Setting either to 0 prevents contouring. (2) Setting either to [] enables contours with automatic selection of levels. (3) Setting either to a positive integer provides a suggestion for the number of levels, with the actual number being set by pretty), which is provided with the integer. (4) Setting either to a vector of numbers specifies those numbers as the levels.

By default, a freezing-point line is drawn (if it is within the range of the data). This behaviour is skipped if draw_freezing is false.

By default, axis names are written in long form; set abbreviate=true for shorter versions.

Information about the analysis is printed if debug is set to true.

Apart from that, the other parameters have the usual meanings for Julia plots. For example, color is set to black, to override the Julia default, etc. In addition to those parameters, the kwargs... argument represents any other argument that is accepted by plot. This is illustrated in the Examples.

Note that specifying seriestype=:line will yield a warning suggesting to use :path instead.

Arguments

  • ctd a Ctd value for which a temperature-salinity diagram will be plotted.

Keywords

  • sigma0_levels a specification of sigma0 values to be contour. If this is an empty vector (which is the default) then the levels are selected automatically by providing pretty with values inferred from ctd. If sigma0_levels equals 0 then no contours are drawn. If it is a positive integer, then it is taken as a suggestion for the number of levels. And, finally, if it is a vector, then it is taken as a specification of the levels to be contoured.

  • sigma0_levels as sigma0_levels, but for spiciness0 contours.

  • draw_freezing a Bool indicating whether to draw a freezing-point curve.

  • abbreviate a Bool indicating whether to abbreviate the axis labels,.

  • debug indicator of debugging level. If this exceeds 0, some information is printed during processing.

  • kwargs... is passed to plot(), to permit further customization; see https://docs.juliaplots.org/stable/ for more information on possibilities.

using OceanAnalysis, Plots, Dates
pkgdir = dirname(dirname(pathof(OceanAnalysis)))
f = joinpath(pkgdir, "data", "ctd.cnv")
ctd = read_ctd_cnv(f);
# Example 1: set title
plot_TS(ctd, title="Built-in CTD file")
# Example 2: just symbols, with no line
plot_TS(ctd, seriestype=:scatter)
# Example 3: just a line, with no symbols
plot_TS(ctd, marker=:none)

See also plot_profile.

source
OceanAnalysis.plot_topographyFunction
plot_topography(topo::Topography;
    xlims=:auto, ylims=:auto, tickdirection=:out,
    domain=:sea, color=:land_sea, clim=:auto,
    draw_coastline=true, land_color=:bisque3, sea_color=:lightblue,
    debug::Int64=0, kwargs...)

Draw a heatmap image of topography.

The domain argument tells whether to display both land and sea values, or just land, or just sea. The default is to plot just the sea, with land a light brown color. The aspect_ratio argument should not be specified as part of kwargs..., because this function sets a reasonable default, based on the latitude at the centre of the plot.

# Waters near Prince Edward Island, Canada
using OceanAnalysis
topo_file = get_topography(-64.8, -61.5, 45.6, 47.2, resolution=1)
topo = read_topography(topo_file)
plot_topography(topo)
source
OceanAnalysis.pressure_from_depthFunction
Compute sea pressure (dbar) from depth (m) and latitude (deg).

Examples

julia> using OceanAnalysis

julia> pressure_from_depth(10.0)
10.082069761243858
source
OceanAnalysis.pressure_from_zFunction
Compute sea pressure (dbar) from vertical coordinate (m) and latitude (deg).

Examples

julia> using OceanAnalysis

julia> pressure_from_z(-10.0)
10.082069761243858
source
OceanAnalysis.prettyFunction
pretty(x, n::Int64=5; debug::Int64=0)

Calculate sub-intervals with 125 scaling

This function finds intervals that are multiples of 1, 2 or 5. The results are useful for contour intervals, because the built-in contour() function (Reference 1) simply divides the range into equal intervals, so that labelled contours can be quite ugly.

Examples

julia> using OceanAnalysis

julia> pretty([22.299, 25.091])
8-element Vector{Float64}:
 22.0
 22.5
 23.0
 23.5
 24.0
 24.5
 25.0
 25.5

References

  1. https://github.com/JuliaGeometry/Contour.jl/blob/daad6eb0b1464dbc7e824bf8384cad54a3b76445/src/Contour.jl#L100)
source
OceanAnalysis.read_adp_rdiFunction
read_adp_rdi(filename::String, ensembles::Union{Int64,StepRange{Int64,Int64},Vector{Int64}}=0; debug::Int64=0)

Read acoustic-Doppler profiler data in RDI "Workhorse-II' format

This function, still in an early phase of development, is designed to read the Teledyne-RDI Workhorse-II PD0 format, as described in Chapter 4 of Reference 1. This format replaces the Workhorse-I PD0 format of Reference 2, which was used as the basis for the read.adp.rdi() function of the R oce package.

It is worth noting that the oce code handles other Teledyne-RDI formats in addition to the Workhorse variety, and it has been tested well with a variety of datasets. It is possible to call R from Julia, so users ought to consider doing so on files that the present function cannot handle.

Arguments

  • filename an ADCP file in the 'PD0' format as described in the Teledyne RD Instruments documentation (Reference 1).

  • ensembles an indication of which ensembles (data profiles) to read. This may be an single integer or a vector of integers. In the first case, if ensembles=0 then the whole file is read, otherwise the stated number of ensembles is read (provided that the file holds that number). In the second case, the value of ensembles dictates the indices of ensembles that are to be read. In both cases, the indices are trimmed to be from 1 to the number of ensembles in the file. The default is to read the whole file. and e.g. ensembles=1:10:101 would read ensemble 1, ensemble 11, and so on, up to ensemble 101.

Keywords

  • debug an integer indicating whether to print information during processing. The default value of 0 means to work quietly, and any larger integer indicates to print some information.

Examples

using OceanAnalysis, Plots

# Load a sample file provided with the package
file = joinpath(dirname(dirname(pathof(OceanAnalysis))), "data", "adp_rdi.000");
adp = read_adp_rdi(file);

# Plot a timeseries of heading
plot(adp["time"], adp["heading"],
    ylab="Heading", label=false, framestyle=:box)

# Plot a heatmap of velocity in the first ensemble
heatmap(adp["velocity"][1, :, :], c=cgrad(:RdBu, rev=true))

# Plot a heatmap of velocity in the first bin vs time and distance
heatmap(adp["time"], adp["distance"], adp["velocity"][:,:,1],
    size=(800,600), ylab="Distance [m]", c=:RdBu)

# List other data in the 'adp' object
keys(adp.data)

# See a particular data item
adp["heading"]

# List the metadata in the 'adp' object
keys(adp.metadata)

# See a particular metadata item
adp["frequency"]

References

  1. Teledyne RD Instruments. “Workhorse II Commands and Output Data Format.” November 2025. P/N 957-6156-00. https://www.teledynemarine.com/en-us/support/SiteAssets/RDI/Manuals%20and%20Guides/Workhorse%20II/WorkHorseCommandsandOutputData_Format.pdf.
  2. Teledyne RD Instruments. “Workhorse Commands and Output Data Format.” 2010.
  3. Teledyne RD Instruments. “Acoustic Doppler Current Profiler Principles of Operation: A Practical Primer.” January 2011. https://www.comm-tec.com/Docs/Manuali/RDI/BBPRIME.pdf.
source
OceanAnalysis.read_amsrFunction
read_amsr(filename::String, field::String="SST"; debug=0)

Reads a measurement stream from an AMSR file.

This returns a value of the Amsr type, with metadata containing the filename along with vectors holding the longitude and latitude of the grid and data holding a matrix of the data element with the indicated name. Using name="?" sidesteps this process, instead returning a vector of strings that may be given as name values.

Arguments

  • filename: a string indicating the location of the local file.

  • field: a string used to identify the data field to be extracted. If field="?" then read_amsr returns a vector of strings containing extractable data. Otherwise, if field names one of those items, then read_amsr returns that dataset.

Keywords

  • debug: An indication of whether to print information during processing. The default value of 0 means to work quietly, and any larger integer indicates to print some information.

Examples

# North Atlantic view, using turbo colour scheme
using OceanAnalysis, Plots
f = "~/data/amsr/RSS_AMSR2_ocean_L3_3day_2025-09-07_v08.2.nc";
d = read_amsr(f, "SST");
longitude = d.metadata["longitude"];
latitude = d.metadata["latitude"];
SST = d.data;
heatmap(longitude, latitude, SST, framestyle=:box,
    xlims=(290, 360), ylims=(20, 60),
    aspect_ratio=1/cos(pi*40/180),
    color=:turbo, size=(800, 550), dpi=300,
    title=f * ": SST", titlefontsize=9)
cl = coastline(:global_fine);
plot!(cl.data.longitude .+ 360, cl.data.latitude,
    seriestype=:shape, color=:bisque3, linewidth=0.8,
    legend=false)
source
OceanAnalysis.read_argoFunction
read_argo(filename::String; column::Int64=1, debug::Int64=0)

Read an Argo file.

Arguments

  • filename a String holding the name of a NetCDF file that holds Argo data.

Keywords

  • column an integer, indicating which profile to read from the file.

  • debug indicator of debugging level. If this exceeds 0, some information is printed during processing.

Return value

The read_argo() function returns an Argo object that has two components, a Dict named .metadata and DataFrame named .data. The .metadata entries are named "cycle", "data_mode", "date_creation", "filename", "latitude", "longitude", "platform", and "time". The .data columns are taken from the source file.

Examples

julia> using OceanAnalysis, Plots

julia> pkgdir = dirname(dirname(pathof(OceanAnalysis)));

julia> f = joinpath(pkgdir, "data", "D4902911_095.nc");

julia> d = read_argo(f);

julia> d.metadata["time"]
2019-10-14T23:43:44.003

julia> d.metadata["latitude"]
40.45216

julia> d.metadata["longitude"]
-66.38298

julia> size(d.data)
(1014, 15)
source
OceanAnalysis.read_argo_indexFunction
read_argo_index(filename::String; trim::Bool=true, header::Int64=9, debug::Int64=0)

Read a file downloaded by get_argo_index.

This relies on there being exactly header lines of header, the last of which names the columns. The default value of 9 works with index files downloaded from the ifremer.fr server, as of 2025-09-08.

First, the date column is converted to a DateTime column named time. If trim is true, then the original date column is removed, along with the the columns named institution, date_update, ocean, and profiler_type.

Return

read_argo_index returns a DataFrame with column names "file", "latitude", "longitude", and "time". Note that the "file" column holds information on the location on remote servers, as is required for use as the file argument of get_argo.

source
OceanAnalysis.read_ctd_cnvFunction
read_ctd_cnv(filename::String; rename::Bool=true, add_teos::Bool=true, debug::Int64=0)

Read a Seabird CTD file in cnv format, optionally adding TEOS-10 variables.

Returns a Ctd object that holds metadata and data. The metadata item is a Dict that holds header (a vector of strings, one per line from the start down to a line containing #END), plus some particular items scanned from that header, e.g. "longitude" and "latitude". The data item is a DataFrame holding the columnar data read from the file. If rename=true, then rename_data is used to rename some of the columns in data to better match oceanographic conventions (e.g. "pr" becomes "pressure"). If the data file indicates temperature is on the T68 scale, then this is converted to the standard modern scale, T90, before saving as temperature.

A message is printed if no data in the file are labelled with names that are recognized as salinity, temperature, or pressure, because these quantities are required for any meaningful CTD dataset.

The value of add_teos is passed to as_ctd, where it indicates whether to add TEOS-10 variables such as SA, CT, sigma0 and spiciness0 to the data portion of the return value.

Examples

julia> using OceanAnalysis

julia> f = joinpath(dirname(dirname(pathof(OceanAnalysis))), "data", "ctd.cnv");

julia> d = read_ctd_cnv(f, add_teos=false);

julia> d.metadata["time"] # note the erroneous year
1903-10-15T11:38:38

julia> d.metadata["latitude"]
44.684266666666666

julia> d.metadata["longitude"]
-63.643883333333335

julia> names(d.data)
8-element Vector{String}:
 "salinity"
 "temperature"
 "pressure"
 "scan"
 "time_seconds"
 "depS"
 "t068"
 "flag"
source
OceanAnalysis.read_ctd_exchangeFunction
read_ctd_exchange(filename::String; add_teos=true, debug::Int64=0)

Read a CTD file in 'exchange' format

Returns a Ctd object that holds metadata (a Dict) and data (a DataFrame). The metadata item is a Dict holding header (the information at the start of the file, down to a line starting DBAR,), along with various quantities, renamed to lower case and some computed quantities, e.g. datetime is constructed from DATE and TIME in the file. The data item is a DataFrame holding the data. The names are altered to be easier to guess, e.g. CTDPRS becomes pressure. Since metadata["header"] holds the original header, users ought to be able to deal with any confusion about names without too much difficulty.

The value of add_teos is passed to as_ctd, where it indicates whether to add TEOS-10 variables such as SA, CT, sigma0 and spiciness0 to the data portion of the return value.

NOTE: The 'exchange' format is commonly supplied by data repositories such as https://cchdo.ucsd.edu, with web links that are typically accompanied with a description containing the word "exchange". The present function cannot handle a related file format, called 'WOCE', but this limitation should not be too problematic because servers tend to supply both formats. (The function reports an error if it is provided with a file in 'WOCE' format.)

Examples

julia> f = joinpath(dirname(dirname(pathof(OceanAnalysis))), "data", "ar07_74JC20140606_00234_00001_ct1.csv");

julia> d = read_ctd_exchange(f);

julia> println(keys(d.metadata))
["latitude", "time", "header", "section", "longitude", "bottom_depth", "station", "expocode", "cast"]

julia> println(first(d.data, 3))
3×8 DataFrame
 Row │ pressure  pressure_flag  temperature  temperature_flag  salinity  salinity_flag  oxygen   oxygen_flag
     │ Float64   Int64          Float64      Int64             Float64   Int64          Float64  Int64
─────┼───────────────────────────────────────────────────────────────────────────────────────────────────────
   1 │      5.0              2      13.1134                 2   34.1855              2    266.0            2
   2 │      7.0              2      13.1195                 2   34.1844              2    266.8            2
   3 │      9.0              2      13.1152                 2   34.1837              2    266.4            2

```

source
OceanAnalysis.read_ctd_rskFunction
read_ctd_rsk(filename::String; add_teos::Bool=true,
    atmospheric_pressure=missing, longitude::Real=-60.0, latitude::Real=40.0,
    debug::Int64=0)

Read a CTD file from an RBR instrument.

These files store data in sqlite format, which the present function handles with the SQLite package. Only a subset of the tables within RBR files are read in this version of the function. The main such tables are: data, which holds the data; channels, which is used to rename the elements in data; and (if present) geodata, which may hold information on the sampling location. Other tables are consulted to learn things like the serial number of the instrument. Some information on the reading process is printed if you call the function with debug=1.

Note that RBR files typically record 'gauge' pressure, which is the sum of atmospheric pressure and the sea pressure. Since oceanographic calculations are typically formulated in terms of sea pressure, so the pressure data stored in RBR files is not inserted into the value returned by this function. Instead, the returned value is set to pressure stored in the file minus atmospheric pressure. The value of atmospheric pressure is taken from the data file if it holds a table called deriveDepth. If such a table is not found, and if the atmosphericPressure argument was not supplied, a default of 10.1325 dbar is used. However, if that argument is supplied, then it supercedes a value stored in the file,

If add_teos is true, then the TEOS10 quantities CT, SA, sigma0 and spiciness0 are computed. These require knowledge of the sampling location, which is inferred from the geodata table (if it exists) or from the function arguments named longitude and latitude, otherwise. (The default values of these parameters specify a position in the western North Atlantic.)

Note that these data files are in a raw form from the instrument, and so some trimming to the downcast portion(s) may be needed. Such features are not yet provided in this package.

Examples

using OceanAnalysis, Plots
ctd = read_ctd_rsk("~/git/oce/create_data/rsk/060130_20150904_1159.rsk");
Sp = plot_profile(ctd, which="salinity");
Tp = plot_profile(ctd, which="temperature");
TS = plot_TS(ctd);
pt = plot(ctd.data.time, ctd.data.pressure);
plot(Sp, Tp, pt, TS, layout=(2,2))
source
OceanAnalysis.read_echosounderFunction
read_echosounder(filename::String; channel::Int64=1, beam::Symbol=:single_beam, sound_speed::Float64=1490.3, debug::Int=0)

Read data from a Biosonics scientific echosounder.

This function is still in development. The read.echosounder() function of the oce R package is being used as a guide, along with the Biosonics document (Reference 1) that describes the file format. FIXME: when code is stabilized, write more here.

Arguments

  • filename string naming the file to be read. It must be in Biosonics DT4 format (reference 1).

Keywords

  • channel an Int64 giving the channel number to read. The default is 1. In the file named in the Examples section, there are two channels, numbered 1 and 2.

  • beam a Symbol indicating which type of beam is sought. In the present version, only :single_beam is permitted. In a later version, :dual_beam and :split_beam may also be permitted.

  • sound_speed a Float64 value indicating the sound speed. The default is to use the sound speed at practical salinity 35, in-situ temperature 10 °C and pressure 30 dbar. (Sound speed is also stored in the files, but Reference 1 says that they ought not to be used.)

  • debug an Int64 value indicating whether to print messages during processing. By default, this is 0, meaning to work quietly.

Examples

# This test will only work if a particular file exists
using OceanAnalysis
filename = "/Users/kelley/Dropbox/data/archive/sleiwex/2008/fielddata/2008-07-01/Merlu/Biosonics/20080701_163942.dt4"
if isfile(filename)
    e = read_echosounder(filename);
end

References

1.BioSonics Advanced Digital Hydroacoustics. “DT4 Data File Format Specification.” BioSonics, May 2017. This is available (after registratration) online at https://www.biosonicsinc.com/support/customer-downloads/

source
OceanAnalysis.read_nonnaFunction
read_nonna(filename::String)

Read a NONNA (NON-NAvigational) bathymetric file.

The Canadian Hydrographic Service provides access to high-resolution bathymetric data for some Canadian waters at https://data.chs-shc.ca/dashboard/map. (This is a GUI-oriented site, and it is somewhat challenging to navigate.) Files are available at both 10-m and 100-m resolution and in a variety of formats. The present function handles only the GeoTIFF format.

Return Value

This returns a Nonna object that holds metadata and data. The metadata item is a Dict that holds easting, northing (both in metres) and some other elements. The data item is a Matrix of the height (in metres) above a nominal sea-level surface.

Examples

using OceanAnalysis, Plots
# Region near East Lawrencetown Beach Provincial Park
n = read_nonna(expanduser("~/data/nonna/NONNA10_4460N06340W.tiff"));
heatmap(n["longitude"], n["latitude"], n.data, color=:turbo,
    framestyle=:box, tickdirection=:out, dpi=300,
    aspect_ratio=1.0, size=(800, 800))
source
OceanAnalysis.read_sectionFunction
read_section(Section::s, debug::Int64=0)

Read an oceanographic section, as downloaded with get_section.

This uses read_ctd_exchange to read WOCE-format CTD files in the directory named dir. These are then aggregated into a section using as_section.

Arguments

  • dir: String naming a directory that holds WOCE-format CTD files.

Keywords

  • debug: an optional value that, if it exceeds 0, indicates that debugging output should be printed during processing.

Examples

using OceanAnalysis, Plots
url = "https://cchdo.ucsd.edu/data/11852/ar07_74JC20140606_ct1.zip"
dir = get_section(url)
section = read_section(dir);
display(section.metadata)
println("Section contains ", length(section.data), " CTD profiles")
source
OceanAnalysis.read_topographyFunction
read_topography(filename::String; debug::Int64 = 0)

Read a topography file that is in NetCDF format. The return value stores longitude in rval.metadata["longitude"], latitude in rval.metadata["latitude"], and depth in rval.data. Note that the depth matrix is transposed, to make it easier to plot.

See also get_topography.

Examples

# Plot world view of ocean depth
using OceanAnalysis, Plots
topo_file = get_topography(:global_coarse);
topo = read_topography(topo_file);
water_depth = -topo.data / 1000.0; # depth (i.e. negative height) in km
water_depth[water_depth .< 0.0] .= NaN; # trim land
heatmap(topo.metadata["longitude"], topo.metadata["latitude"], water_depth,
    asp=1.0, framestyle=:box, xlims=[-180,180], ylims=[-90,90],
    color=cgrad(:deep, rev=false), dpi=300)
cl = coastline();
plot!(cl.data.longitude, cl.data.latitude, color=:black, legend=false, linewidth=0.5)
source
OceanAnalysis.rename_dataFunction
rename_data(names)

Rename data items from labels used in files to names used in code.

Arguments

  • names a String holding a name to be converted, or a vector of such strings.

  • number_replicates a Bool value indicating whether to prevent duplicated names by appending numbers to duplicates. This is true by default. For example, the first salinity would be named salinity, while the second would be named salinity2.

Return value

A String or vector of String items, holding new names. If any of the converted names appear more than once, then digits are appended (see last example).

Examples

julia> using OceanAnalysis
julia> rename_data("CTDPRS")
"pressure"

julia> rename_data(["CTDPRS", "CTDTMP"])
2-element Vector{String}:
 "pressure"
 "temperature"

julia> rename_data(["CTDPRS", "CTDTMP", "CTDTMP_FLAG"])
3-element Vector{String}:
 "pressure"
 "temperature"
 "temperature_flag"
source
OceanAnalysis.running_meanFunction
running_mean(x::vector{float64}, k::int64=3)

Compute running mean of a vector x, over a window of width k.

Examples

i = 1:100
x = sin.(2 * pi * i / 50)
ymean = running_mean(x, 3)
scatter(x, label="data", ms=2)
plot!(ymean, label="mean")
source
OceanAnalysis.running_medianFunction
running_median(x::vector{float64}, k::int64=3)

Compute running median of a vector x, over a window of width k.

Examples

i = 1:100
x = sin.(2 * pi * i / 50)
ymedian = running_median(x, 3)
scatter(x, label="data", ms=2)
plot!(ymedian, label="median")
source
OceanAnalysis.SAFunction
SA(salinity, pressure, longitude, latitude)
SA(ctd)

Compute Absolute Salinity (SA), using gsw_sa_from_sp() in the GibbsSeaWater package.

The first form takes single values and returns a single value. The second form extracts values from a Ctd object and then calls the first form as SA.() so that it returns a vector of SA values.

Examples

julia> using OceanAnalysis

julia> SA(35.0, 100.0, -30.0, 30.0)
35.165308620244
source
OceanAnalysis.scale_barFunction
scale_bar(distance::Real=100.0, x=:left, y=:top;
    linewidth::Real=2.0, fontsize::Real=8)

Add a horizontal scalebar to a plot made with [plot_coastline]@ref).

The length of the scalebar, in km, is given by distance, at a position dictated by x and y. The value of x must be :left, :right or a number (for longitude), and the value of y must be :bottom, :top or a number (for latitude). The default is to place the scale bar at the top-left. If none of the corners are suitable, e.g. if the label covers important parts of the plot, use numeric values for x and y as the longitude and latitude of the beginning of the line indicating the scale. The width of the line is given by linewidth, and the fontsize of the label is given by fontsize.

Examples

using OceanAnalysis, Plots
cl = coastline();
plot_coastline(cl, xlims=(-70, -60), ylims=(42, 48))
scale_bar(100.0)
source
OceanAnalysis.set_teosFunction
set_teos(x::OA; debug::Int64=0)

Add, or modify, TEOS-10 components to hydrographic data.

Compute the TEOS-10 quantities SA (Absolute Salinity), CT (Conservative Temperature), sigma0 (potential density anomaly with respect to surface pressure), and spiciness0 (seawater spiciness with respect to surface pressure). These items are inserted into the data component of the returned value. If they are already present in x, then new values are inserted in the return value.

An error is reported if the x.data lacks salinity, temperature or pressure, or if x.metadata lacks longitude or latitude.

source
OceanAnalysis.section_is_griddedFunction
section_is_gridded(Section:section; debug:Int64=0)

Return true if section is gridded (that is, if it has more than 1 CTD station, and if the pressure levels match across all the CTD stations.

source
OceanAnalysis.station_mapFunction
station_map(longitude, latitude; scale::Real=5.0, debug::Int64=0, kwargs...)

Using plot_coastline, draw a map that shows the location of a station (or stations) specified by longitude and latitude, each of which may be a single number or a vector of numbers (with longitude in the -180 to +180 convention). The map span is computed automatically by computing the distance between the centroid of the stations and the nearest point of land and also computing the span across the stations. The maximum of these two distances is multiplied by scale, and from this the x and y limits of the plot are set. Altering the value of scale is thus the way a user can control the view. The station locations are drawn by calling scatter from the Plots package, to which the kwargs... elements are passed directly; the example shows how to use this fact to alter the station symbols.

Map projections are not offered by station_map; to get such views, consider using the GMT package.

Examples

using OceanAnalysis, Plots
# Red circle marks a station south of due east of Fort Louisbourg
# and due south of Saint Pierre and Miquelon.
p1 = station_map(-56.33, 45.90)
# The same, but using a large light-blue diamond
p2 = station_map(-56.33, 45.90,
    markercolor=:lightblue, markershape=:diamond, markersize=6)
plot(p1, p2)
source
OceanAnalysis.subset_amsrFunction
subset_amsr(a::Amsr, lonlims, latlims; debug::Int64=0)

Subset an Amsr object to a specified longitude and latitude range.

Arguments

  • a: an Amsr object.

  • lonlims: A numeric tuple of length 2 specifying the minimum and maximum longitude values to be retained.

  • latlims: A numeric tuple of length 2 specifying the minimum and maximum latitude values to be retained.

  • debug: An indication of whether to print information during processing. The default value of 0 means to work quietly, and any larger integer indicates to print some information.

source
OceanAnalysis.subset_ctdFunction
subset_ctd(ctd::Ctd, keep_levels::Union{BitVector,Vector{Bool}}; debug::Int64=0)

Subset a CTD object by levels.

This returns a copy of Ctd that has the same metadata, but for which the data holds only rows specified by the logical vector keep_levels.

Examples

using OceanAnalysis, Plots
f = joinpath(dirname(dirname(pathof(OceanAnalysis))), "data",
"D4902911_095.nc")
argo = read_argo(f)
ctd = as_ctd(argo)
a = plot_TS(ctd, title="Original")
ctd_top = subset_ctd(ctd, ctd["pressure"] .< 300)
b = plot_TS(ctd)
plot(a, b, title="Top 300m")
source
OceanAnalysis.subset_ctd!Function
subset_ctd!(ctd::Ctd, keep_levels::Union{BitVector,Vector{Bool}}; debug::Int64=0)

This works in the same way as subset_ctd(<Ctd>), except that the original Ctd is altered in-place.

Examples

using OceanAnalysis, Plots
f = joinpath(dirname(dirname(pathof(OceanAnalysis))), "data",
"D4902911_095.nc")
argo = read_argo(f)
ctd = as_ctd(argo)
a = plot_TS(ctd, title="Original")
subset_ctd!(ctd, ctd["pressure"] .< 300)
b = plot_TS(ctd, title="Original (altered)")
plot(a, b)
source
OceanAnalysis.T90_from_T48Function
T90 = T90_from_T48(T48::Float64)

Convert a temperature from the T48 scale to the T90 scale.

See also T90_from_T68.

Examples

julia> using OceanAnalysis

julia> T90_from_T48(10.0)
9.993641526033752
source
OceanAnalysis.T90_from_T68Function
T90 = T90_from_T68(T68::Float64)

Convert a temperature from the T68 scale to the T90 scale.

See also T90_from_T48.

Examples

julia> using OceanAnalysis

julia> T90_from_T68(10.0)
9.997600575861792
source
OceanAnalysis.tocFunction
toc(x::OA)

Print a table of contents for an OA object.

Examples

julia> using OceanAnalysis

julia> f = joinpath(dirname(dirname(pathof(OceanAnalysis))), "data", "ctd.cnv");

julia> d = read_ctd_cnv(f);

julia> d = read_ctd_cnv(f, add_teos=false);

julia> toc(d)
Ctd object
    metadata: Dict with keys ["latitude", "time", "header", "longitude", "filename"]
    data: DataFrame with columns ["salinity", "temperature", "pressure", "scan", "timeS", "pr", "depS", "t090", "sal00", "flag"]
source
OceanAnalysis.xyz_to_enuFunction
xyz_to_enu(adp::Adp; declination::Float64=0.0, debug::Int64=0)

Change velocity in an RDI Workhorse Adp object from xyz to enu coordinates

This is done by using the heading, pitch and roll vectors that are stored within adp. Note that declination is added to heading, to allow for compass correction. See read_adp_rdi for how to read an RDI Adp object, and beam_to_xyz for how to convert it from beam to xyz coordinates.

Examples

```juliadoc using OceanAnalysis, Plots file = joinpath(dirname(dirname(pathof(OceanAnalysis))), "data", "adprdi.000") beam = readadprdi(file); xyz = beamtoxyz(beam); enu = xyzto_enu(xyz); v = enu["velocity"];

source