OceanAnalysis Functions

Data types:

OceanAnalysis.OAType
Base type in the OceanAnalysis package.

This is an abstract type. The other types in the package will derive from this. At the moment, the only such case is Ctd.

source
OceanAnalysis.AmsrType
A type to hold AMSR data (SUBJECT TO CHANGE)

This holds AMSR satellite data as read by read_amsr.

The metadata element is a Dict that holds the source filename, the longitude and latitude, and the name of the stored variable.

The data element stores the data.

source
OceanAnalysis.CtdType
A type to hold CTD data

Ctd is a type used to store data from CTD instruments and Argo floats. It takes the form of a Struct that derives from the base OA type, and holds 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 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 and read_argo and can be passed to plotting functions plot_profile and plot_TS, and by several 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.

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.

source
OceanAnalysis.TopographyType
A type to hold topography data (SUBJECT TO CHANGE)

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.

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(salinity::Vector{Float64}, temperature::Vector{Float64}, pressure::Vector{Float64};
    longitude::Float64=-63.0, latitude::Float64=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

julia> using OceanAnalysis

julia> as_ctd([32.], [15.], [0.], add_teos=false)
Ctd(Dict{String, Any}("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=-63., latitude=40.)
Ctd(Dict{String, Any}("latitude" => 40.0, "time" => nothing, "longitude" => -63.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.1516  15.0642  23.6653   0.0686905)


julia> as_ctd([32.], [15.], [0.], longitude=-63., latitude=30.)
Ctd(Dict{String, Any}("latitude" => 30.0, "time" => nothing, "longitude" => -63.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.1511  15.0642  23.6649   0.0683062)
source
OceanAnalysis.coastlineFunction
coastline(symbol::Symbol=:world_fine)

Access a built-in coastline dataset. The only valid choices for name are :world_coarse and :world_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)

Return a coastline stored in the named CSV file (in either text or gzipped form).

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::Vector{Real}, latitude::Vector{Real})

Create a Coastline from longitude and latitude values. Use NaN values for each of these 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.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.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_amsr_fileFunction
get_amsr_file(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_file 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_file.

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_file 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_argo_fileFunction
get_argo_file(file::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.

  • 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.

Returns

  • get_argo_file 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_file(index.file[end])
argo = read_argo(argo_file)
plot_profile(argo, which="CT")
source
OceanAnalysis.get_fileFunction
Download a remote file or identify an existing version if it is young

The contents of url are downloaded and stored as file, but only if either file does not exist locally or its age is less than age days.

source
OceanAnalysis.get_nc_valueFunction
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.

  • require_value boolean value indicating whether to report an error if the desired element consists entirely of bad values.

source
OceanAnalysis.get_topography_fileFunction
get_topography_file(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_file(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_file 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.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 bc="nearest" to control what happens near boundaries. The user has no control over these things, although this might change in a future version of N2().

The user's control rests in s, a smoothing parameter that is passed to Dierckx:Spline1D(). If not specified by the user, this defaults to a value that yields N² curves that are similar to those computed with a default call to swN2() in the oce R package. Users may elect to use larger s values for smoother curves, or smaller ones to get more detail. It would be a mistake not to pair experiments with s values with plots. As a start, it might be useful to examine Reference 2, which compares the R and Julia results.

References

  1. https://github.com/JuliaMath/Dierckx.jl
  2. https://github.com/dankelley/OceanAnalysis.jl/issues/13

Examples

# Demonstrate N2()
using OceanAnalysis, Plots
pkgdir = dirname(dirname(pathof(OceanAnalysis)))
filename = joinpath(pkgdir, "data", "ctd.cnv")
ctd = read_ctd_cnv(filename);
# Basic plot (left), profile-plot in oceanographic convention (right)
p1=plot(N2(ctd), z_from_pressure.(ctd.data.pressure))
p2=plot_profile(ctd, which="N2")
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. The first is a Dict named metadata, which has contents that vary with data type. The second is named data, the contents of which depend on the data type. The elements of both may be accessed directly using the dot notation, with for example ctd.data.salinity refering to the salinity column of ctd data read by read_ctd_cnv. It is also possible to retrieve (or set) that value using ctd["salinity"]. The advantage of this notation is that it can locate information whether it is in the data or the metadata component of the object. It is also possible to obtain derived information, e.g. ctd["SA"] calculates and then returns the Absolute Salinity, which is not typically stored in CTD files, but which can be computed from the stored information. Other derivable items include "CT" (Conservative Temperature), "sigma0" (potential density anomaly with respect to surface pressure), and spiciness0 (spiciness).

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;
    xlims=(-180., 180.), ylims=(-90., 90.),
    seriestype=:shape, color=:bisque3, linewidth=0.5, tickdirection=:out,
    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

  • coastline the coastline, as constructed using coastline or (less commonly) Coastline.

  • xlims and ylims control the ranges of the longitude and latitude axes, respectively.

  • seriestype, color and linewidth control the rendering of land regions. These values are passed to the base-level plot function; for details, see the documentation provided by the Plots package.

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

  • coastline the coastline, as constructed using coastline or, by more advanced users, using Coastline.

Keywords

  • seriestype, color and linewidth control the rendering of land regions. These values are passed to the base-level plot function; for details, see the documentation provided by the Plots package.
source
OceanAnalysis.plot_profileFunction
plot_profile(ctd::Ctd; which::String="CT", vertical::String="pressure",
    abbreviate::Bool=false, legend::Bool=false, tickfontsize=8, tickdirection=:out,
    labelfontsize=8, debug::Int64=0, kwargs...)

Plot an oceanographic profile for data contained in ctd, showing how the variable named by which depends on pressure. The variable is drawn on the x axis and pressure on the y axis. Following oceanographic convention, pressure increases downwards on the page and the "x" axis is drawn at the top. The permitted values of which are "CT" for the Gibbs Seawater formulation of Conservative Temperature, "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.

The default Julia font sizes on axes are overridden in this function, with 8-point being used for both the numbers on axes (tickfontize) and the names of axes (labelfontsize). (The tickfontsize matches the Julia default, but the labelfontsize is smaller than the Julia default. The idea is to not waste space by using fonts that are larger than what journals require.)

The kwargs... argument is used for arguments to be sent to plot(). For example, the default way to display the profile diagram is constructed with a blue line connecting points, but using e.g.

plot_profile(ctd, which="SA", seriestype=:scatter, seriescolor=:red)

yields red-filled circles, instead; see https://docs.juliaplots.org/stable/ for more on the many plotting controls available in Julia.

See also the plot_TS function.

Examples

using OceanAnalysis, Plots
# Read an Argo file
pkgdir = dirname(dirname(pathof(OceanAnalysis)))
f = joinpath(pkgdir, "data", "D4902911_095.nc")
d = read_argo(f, 1);
# Plot profiles of Conservative Temperature, Absolute Salinity, and potential
# density anomaly with respect to surface pressure.
p1 = plot_profile(d, which="CT")
p2 = plot_profile(d, which="SA")
p3 = plot_profile(d, which="sigma0")
plot(p1, p2, p3, layout=(1, 3), size=(800, 400))
source
OceanAnalysis.plot_TSFunction
plot_TS(ctd::Ctd; sigma0_levels=[], spiciness0_levels=0,
    draw_freezing=true, abbreviate=false,
    framestyle=:box, color=:black, seriestype=:scatter, markersize=2.0, linewidth=1.0,
    legend=false, gridstyle=:dash, tickfontsize=8, tickdirection=:out,
    labelfontsize=8, 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 example, which a title is added to the plot for a built-in CNV-formatted CTD file.

using OceanAnalysis, Plots, Dates
pkgdir = dirname(dirname(pathof(OceanAnalysis)))
f = joinpath(pkgdir, "data", "ctd.cnv")
ctd = read_ctd_cnv(f);
plot_TS(ctd, title="Built-in CTD file", titlefontsize=9)

See also plot_profile.

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_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, add_teos::Bool=true,
    require_valid::Bool=true, debug::Int64=0)

Read an Argo file and return a Ctd object that holds salinity, temperature, pressure (and derived columns) but no other columns from the file. This function is in an early stage of development; please report problems as issues on <www.github.com/dankelley/OceanAnalysis.jl/issues>.

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.

If require_valid is true (the default) then an error is reported if the file lacks one of three required data columns, or either longitude or latitude. An error is also reported if any of these items consists entirely of missing values. This is because such files are unlikely to be of much use. In some cases, setting require_valid to false may permit the file to be read, but this has not been tested, since the results in such cases are not likely to be of practical use.

Set debug to a positive integer to cause read_argo() to print messages during processing. This can be handy if problems arise.

Return value

The read_argo() function returns a Ctd 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 named "pressure", "salinity" and "temperature", as copied from fields in the NetCDF file named "PRES", "PSAL" and "TEMP"; no other NetCDF fields are copied in this version of read_argo().

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> first(d.data,3)
3×7 DataFrame
 Row │ salinity  temperature  pressure  SA       CT       sigma0   spiciness0
     │ Float64   Float64      Float64   Float64  Float64  Float64  Float64
─────┼────────────────────────────────────────────────────────────────────────
   1 │   34.913       19.513      0.48  35.0786  19.5079  24.8272     3.31464
   2 │   34.91        19.527      1.0   35.0756  19.5219  24.8213     3.31603
   3 │   34.912       19.524      2.0   35.0776  19.5187  24.8237     3.31669
source
OceanAnalysis.read_argo_indexFunction
read_argo_index(file::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_file.

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

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

Returns a Ctd object that holds metadata (a Dict) and data (a DataFrame). metadata item 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. data holds the columnar data read from the file, along with renamed values in standard nomenclature. At present, the only renamed items are salinity, temperature, and pressure. Note that 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"]
2003-10-15T11:38:38

julia> d.metadata["latitude"]
44.684266666666666

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

julia> names(d.data)
10-element Vector{String}:
 "salinity"
 "temperature"
 "pressure"
 "scan"
 "timeS"
 "pr"
 "depS"
 "t090"
 "sal00"
 "flag"
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_file.

Examples

# Plot world view of ocean depth
using OceanAnalysis, Plots
topo_file = get_topography_file(: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.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.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