Skip to content

amundsen

Amundsen Science.

This module contains the parser for the Amundsen INT format. This format is used to store oceanographic data generated by the Amundsen Science and ArcticNet programs.

csv_format(path, encoding='UTF-8', map_to_vocabulary=True, generate_depth=True, separator=',', encoding_error='strict')

Parse Amundsen CSV format.

Parameters:

Name Type Description Default
path str

file path to parse.

required
encoding str

File encoding. Defaults to "Windows-1252".

'UTF-8'
map_to_vocabulary bool

Rename variables to vocabulary. Defaults to True.

True
generate_depth bool

Generate depth variable. Defaults to True.

True
separator str

Separator for the data. Defaults to r",".

','
encoding_error str

Encoding error handling. Defaults to "strict".

'strict'

Returns:

Type Description
Dataset

xr.Dataset

Source code in ocean_data_parser/parsers/amundsen.py
def csv_format(
    path: str,
    encoding: str = "UTF-8",
    map_to_vocabulary: bool = True,
    generate_depth: bool = True,
    separator: str = ",",
    encoding_error="strict",
) -> xr.Dataset:
    """Parse Amundsen CSV format.

    Args:
        path (str): file path to parse.
        encoding (str, optional): File encoding. Defaults to "Windows-1252".
        map_to_vocabulary (bool, optional): Rename variables to vocabulary. Defaults to True.
        generate_depth (bool, optional): Generate depth variable. Defaults to True.
        separator (str, optional): Separator for the data. Defaults to r",".
        encoding_error (str, optional): Encoding error handling. Defaults to "strict".

    Returns:
        xr.Dataset
    """
    return int_format(
        path=path,
        encoding=encoding,
        map_to_vocabulary=map_to_vocabulary,
        generate_depth=generate_depth,
        separator=separator,
        encoding_error=encoding_error,
    )

int_format(path, encoding='Windows-1252', map_to_vocabulary=True, generate_depth=True, separator='\\s+', encoding_error='strict')

Parse Amundsen INT format.

Parameters:

Name Type Description Default
path str

file path to parse.

required
encoding str

File encoding. Defaults to "Windows-1252".

'Windows-1252'
map_to_vocabulary bool

Rename variables to vocabulary. Defaults to True.

True
generate_depth bool

Generate depth variable. Defaults to True.

True
separator str

Separator for the data. Defaults to r"\s+".

'\\s+'
encoding_error str

Encoding error handling. Defaults to "strict".

'strict'

Returns:

Type Description
Dataset

xr.Dataset

Source code in ocean_data_parser/parsers/amundsen.py
def int_format(
    path: str,
    encoding: str = "Windows-1252",
    map_to_vocabulary: bool = True,
    generate_depth: bool = True,
    separator: str = r"\s+",
    encoding_error="strict",
) -> xr.Dataset:
    r"""Parse Amundsen INT format.

    Args:
        path (str): file path to parse.
        encoding (str, optional): File encoding. Defaults to "Windows-1252".
        map_to_vocabulary (bool, optional): Rename variables to vocabulary. Defaults to True.
        generate_depth (bool, optional): Generate depth variable. Defaults to True.
        separator (str, optional): Separator for the data. Defaults to r"\s+".
        encoding_error (str, optional): Encoding error handling. Defaults to "strict".

    Returns:
        xr.Dataset
    """
    metadata = default_global_attributes.copy()

    # Ignore info.int files
    if path.endswith("_info.int"):
        logger.warning("Ignore *_info.int files: {}", path)
        return

    logger.debug("Read {}", path)
    with open(path, encoding=encoding, errors=encoding_error) as file:
        # Parse header
        for header_line_idx, line in enumerate(file):
            line = line.rstrip()
            if re.match(r"^%\s*$", line) or not line:
                continue
            elif (
                line and not re.match(r"\s*%", line) and (line[0] == " " or "," in line)
            ):
                break
            elif ":" in line:
                key, value = line.strip()[1:].split(":", 1)
                metadata[key.strip()] = value.strip()
            elif line in VARIABLE_MAPPING_FIXES:
                metadata[VARIABLE_MAPPING_FIXES[line]] = line[2:]
            elif re.match(r"% .* \[.+\]", line):
                logger.warning(
                    "Unknown variable name will be saved to unknown_variables_information: {}",
                    line,
                )
                metadata["unknown_variables_information"] += line + "\n"
            else:
                logger.warning("Unknown line format: {}", line)

        # Review metadata
        if metadata == default_global_attributes:
            logger.warning("No metadata was captured in the header of the INT file.")

    # Parse data
    if separator == r"\s+":
        # Fix problematic variable names
        for name, fixed_name in REPLACE_VARIABLES.items():
            line = line.replace(name, fixed_name)
        names = re.split(r"\s+", line.strip())
    elif separator == ",":
        names = line.strip().split(",")

    if len(set(names)) != len(names):
        duplicated = {name: n for name, n in Counter(names).items() if n > 1}
        logger.warning("Duplicated variable names detected: {}", duplicated)
        # Add index to duplicate names
        new_names = []
        for name in names:
            if name not in new_names:
                new_names.append(name)
            else:
                new_names.append(f"{name}_{new_names.count(name)}")
        names = new_names

    df = pd.read_csv(
        path,
        encoding=encoding,
        header=header_line_idx,
        skiprows=[header_line_idx] if separator == r"\s+" else [],
        sep=separator,
        names=names,
        encoding_errors=encoding_error,
    )
    if len(df.columns) != len(names):
        raise ValueError(
            f"Number of columns ({len(df.columns)}) doesn't match the number of variables ({len(names)})"
        )

    # Sort column attributes
    variables = _extract_variable_attributes_from_header(metadata, df.columns)
    if "Date" in df and "Hour" in df:
        df = _convert_timestamp(df)

    # Convert to xarray object
    ds = df.to_xarray()

    # Standardize global attributes
    ds.attrs = {
        _standardize_attribute_name(name): _standardize_attribute_value(
            value, name=name
        )
        for name, value in metadata.items()
    }

    # Fix global attribute
    ds = _fix_station_global_attribute(ds)

    # Generate instrument_depth variable
    pressure = [var for var in ds if var in ("Pres", "PRES")]
    if (
        generate_depth
        and pressure
        and ("Lat" in ds or "initial_latitude_deg" in ds.attrs)
    ):
        logger.info(
            "Generate instrument_depth from TEOS-10: -1 * gsw.z_from_p(ds['Pres'], {})",
            "ds['Lat']" if "Lat" in ds else "ds.attrs['initial_latitude_deg']",
        )
        ds["instrument_depth"] = -z_from_p(
            ds[pressure[0]],
            ds["Lat"] if "Lat" in ds else ds.attrs["initial_latitude_deg"],
        )

    # Map variables to vocabulary
    file_type = _get_file_type(path)
    variables_vocabulary = amundsen_vocabulary(file_type)
    variables_to_rename = {}
    for var in ds:
        if var not in variables or re.sub(r"_\d+$", "", var) in IGNORED_VARIABLES:
            continue

        ds[var].attrs = variables[var]
        if "long_name" in ds[var].attrs:
            ds[var].attrs["long_name"] = ds[var].attrs["long_name"].strip()

        # Include variable attributes from the vocabulary
        if not map_to_vocabulary:
            continue
        elif var not in variables_vocabulary:
            logger.warning(
                "No vocabulary is available for variable mapping '{}': {}",
                var,
                ds[var].attrs,
            )
            continue

        # Match vocabulary
        var_units = ds[var].attrs.get("units")
        for item in variables_vocabulary[var]:
            accepted_units = item.get("accepted_units")
            if (
                var_units is None  # Consider first if no units
                or var_units == item.get("units")
                or (
                    accepted_units
                    and (
                        re.fullmatch(accepted_units, var_units)
                        or "None" in accepted_units.split("|")
                    )
                )
            ):
                if "rename" in item:
                    variables_to_rename[var] = item["rename"]

                ds[var].attrs = {
                    key: value
                    for key, value in item.items()
                    if key not in ["accepted_units", "rename", "file_type"]
                }
                if "apply_function" in ds[var].attrs:
                    ds = apply_function(ds, var)
                break
        else:
            logger.warning(
                "No Vocabulary available for {} [{}]: {}",
                var,
                var_units,
                str(ds[var].attrs),
            )

    # Review rename variables
    already_existing_variables = {
        var: rename for var, rename in variables_to_rename.items() if rename in ds
    }
    if already_existing_variables:
        logger.error(
            "Can't rename variable {} since it already exist",
            already_existing_variables,
        )

    if variables_to_rename:
        logger.info("Rename variables: {}", variables_to_rename)
        ds = ds.rename(variables_to_rename)

    # Assign dimensions only after renaming variables
    ds = _assign_dimensions(ds, instrument=file_type)

    # Standardize dataset to be compatible with ERDDAP and NetCDF Classic
    ds = standardize_dataset(ds)
    return ds

Rosette Vocabulary

variable_name long_name units standard_name sdn_parameter_urn sdn_uom_urn alternative_label rename sdn_parameter_name sdn_uom_name accepted_units apply_function Units
Pres Sea Pressure decibars sea_water_pressure_due_to_sea_water SDN:P01::PRESPR01 SDN:P06::UPDB Pres_Z PRES Pressure (spatial coordinate) exerted by the water body by profiling pressure sensor and correction to read zero at sea level Decibars nan nan nan
PresPar nan nan nan nan nan nan nan nan nan nan nan nan
PRES Sea Pressure decibars sea_water_pressure_due_to_sea_water SDN:P01::PRESPR01 SDN:P06::UPDB Pres_Z nan Pressure (spatial coordinate) exerted by the water body by profiling pressure sensor and correction to read zero at sea level Decibars nan nan nan
DEPH Depth m depth SDN:P01::ADEPZZ01 SDN:P06::ULAA Depth depth Depth (spatial coordinate) relative to water surface in the water body Metres nan nan nan
Temp Temperature (ITS-90) degC sea_water_temperature SDN:P01::TEMPS901 SDN:P06::UPAA CTDTmp90 TE90 Temperature (ITS-90) of the water body by CTD or STD Degrees Celsius degC deg C degrees C
TE90 Temperature (ITS-90) degC sea_water_temperature SDN:P01::TEMPS901 SDN:P06::UPAA CTDTmp90 nan Temperature (ITS-90) of the water body by CTD or STD Degrees Celsius degC deg C degrees C
Sal Practical Salinity PSU sea_water_practical_salinity SDN:P01::PSALST01 SDN:P06::PSUX P_sal_CTD PSAL Practical salinity of the water body by CTD and computation using UNESCO 1983 algorithm Practical salinity units PSU psu 1e-3
PSAL Practical Salinity PSU sea_water_practical_salinity SDN:P01::PSALST01 SDN:P06::PSUX P_sal_CTD nan Practical salinity of the water body by CTD and computation using UNESCO 1983 algorithm Practical salinity units PSU psu 1e-3
TRB_ Turbidity NTU sea_water_turbidity SDN:P01::TURBXXXX SDN:P06::USTU Turbidity TRB Turbidity of water in the water body Nephelometric Turbidity Units NTU nan nan
TURB Turbidity FTU sea_water_turbidity SDN:P01::TURBXXXX SDN:P06::UFNU Turbidity nan Turbidity of water in the water body Formazin Nephelometric Units FTU nan nan
DateTime nan YYYY-MM-DDThh🇲🇲ss.sssZ nan SDN:P01::DTUT8601 SDN:P06::UYMD time time Date and time (UT in ISO8601 format to known precision) Gregorian date (yyyymmdd) nan nan nan
Lat Latitude degrees_north latitude SDN:P01::ALATZZ01 SDN:P06::DEGN nan latitude Latitude north Degrees north degrees N nan nan
Long Longitude degrees_east longitude SDN:P01::ALONZZ01 SDN:P06::DEGE Lon longitude Longitude east Degrees east degrees E nan nan
Lon Longitude degrees_east longitude SDN:P01::ALONZZ01 SDN:P06::DEGE Lon longitude Longitude east Degrees east degrees E nan nan
Trans Light Transmission % nan SDN:P01::POPTDR01 SDN:P06::UPCT Trans_Red_25cm TRAN Transmittance (red light wavelength) per 25cm of the water body by 25cm path length red light transmissometer Percent nan nan nan
Tran Light Transmission % nan SDN:P01::POPTDR01 SDN:P06::UPCT Trans_Red_25cm TRAN Transmittance (red light wavelength) per 25cm of the water body by 25cm path length red light transmissometer Percent nan nan nan
TRAN Light Transmission % nan SDN:P01::POPTDR01 SDN:P06::UPCT Trans_Red_25cm nan Transmittance (red light wavelength) per 25cm of the water body by 25cm path length red light transmissometer Percent nan nan nan
Fluo chl-a fluorescence ug/L mass_concentration_of_chlorophyll_a_in_sea_water SDN:P01::CPHLPR01 SDN:P06::UGPL chla_water_ISfluor FLOR Concentration of chlorophyll-a {chl-a CAS 479-61-8} per unit volume of the water body [particulate >unknown phase] by in-situ chlorophyll fluorometer Micrograms per litre ug/L mg/m\^3 nan
FLOR chl-a fluorescence ug/L mass_concentration_of_chlorophyll_a_in_sea_water SDN:P01::CPHLPR01 SDN:P06::UGPL chla_water_ISfluor nan Concentration of chlorophyll-a {chl-a CAS 479-61-8} per unit volume of the water body [particulate >unknown phase] by in-situ chlorophyll fluorometer Micrograms per litre ug/L mg/m\^3 mg/m**3'
DOXY Dissolved Oxygen uM mole_concentration_of_dissolved_molecular_oxygen_in_sea_water SDN:P01::DOXYZZ01 SDN:P06::UPOX WC_dissO2_IS OXYM Concentration of oxygen {O2 CAS 7782-44-7} per unit volume of the water body [dissolved plus reactive particulate phase] by in-situ sensor Micromoles per litre ml/l x*44.661 nan
O2 Dissolved Oxygen uM mole_concentration_of_dissolved_molecular_oxygen_in_sea_water SDN:P01::DOXYZZ01 SDN:P06::UPOX WC_dissO2_IS OXYM Concentration of oxygen {O2 CAS 7782-44-7} per unit volume of the water body [dissolved plus reactive particulate phase] by in-situ sensor Micromoles per litre uM µM µM
O2 Dissolved Oxygen uM mole_concentration_of_dissolved_molecular_oxygen_in_sea_water SDN:P01::DOXYZZ01 SDN:P06::UPOX WC_dissO2_IS OXYM Concentration of oxygen {O2 CAS 7782-44-7} per unit volume of the water body [dissolved plus reactive particulate phase] by in-situ sensor Micromoles per litre ml/l x*44.661 nan
OXYM Dissolved Oxygen uM µM µM mole_concentration_of_dissolved_molecular_oxygen_in_sea_water SDN:P01::DOXYZZ01 SDN:P06::UPOX WC_dissO2_IS nan Concentration of oxygen {O2 CAS 7782-44-7} per unit volume of the water body [dissolved plus reactive particulate phase] by in-situ sensor Micromoles per litre
ASAL Absolute Salinity (TEOS-10) g/kg sea_water_absolute_salinity SDN:P01::ASLTZZ01 SDN:P06::UGKG A_sal nan Absolute salinity of the water body Grams per kilogram nan nan nan
Asal Absolute Salinity TEOS-10 g/kg sea_water_absolute_salinity SDN:P01::ASLTZZ01 SDN:P06::UGKG A_sal ASAL Absolute salinity of the water body Grams per kilogram nan nan nan
Sigt Sigma-T kg/m^3 sea_water_sigma_t SDN:P01::SIGTEQST SDN:P06::UKMC SigmaT SIGT Sigma-T of the water body by computation from salinity and temperature using UNESCO algorithm Kilograms per cubic metre nan nan nan
SIGT Sigma-T kg/m^3 sea_water_sigma_t SDN:P01::SIGTEQST SDN:P06::UKMC SigmaT nan Sigma-T of the water body by computation from salinity and temperature using UNESCO algorithm Kilograms per cubic metre nan nan nan
sigt Sigma-T kg/m^3 sea_water_sigma_t SDN:P01::SIGTEQST SDN:P06::UKMC Sigma-T SIGT Sigma-T of the water body by computation from salinity and temperature using UNESCO algorithm Kilograms per cubic metre nan nan nan
sigthe Sigma-Theta (Rho(S, Theta, 0)-1000) kg/m^3 sea_water_sigma_theta SDN:P01::SIGTEQ01 SDN:P06::UKMC SigmaTheta SIGP Sigma-theta of the water body by computation from salinity and potential temperature using UNESCO algorithm Kilograms per cubic metre nan nan nan
Sigthe Sigma-Theta (Rho(S, Theta, 0)-1000) kg/m^3 sea_water_sigma_theta SDN:P01::SIGTEQ01 SDN:P06::UKMC SigmaTheta SIGP Sigma-theta of the water body by computation from salinity and potential temperature using UNESCO algorithm Kilograms per cubic metre nan nan nan
Dens Sigma-T (Rho (S, T, 0) - 1000) kg/m^3 sea_water_sigma_t nan nan nan nan nan nan kg/m^3 kg/m**3 nan
N2 Brunt-Vaisala Frequency s^-2 square_of_brunt_vaisala_frequency_in_sea_water nan nan nan nan nan nan nan nan nan
Theta Potential Temperature degC sea_water_potential_temperature nan nan theta nan nan nan degC deg C nan
theta Potential Temperature degC sea_water_potential_temperature nan nan nan nan nan nan nan nan nan
CDOM CDOM fluorescence mg/m^3 concentration_of_colored_dissolved_organic_matter_in_sea_water_expressed_as_equivalent_mass_fraction_of_quinine_sulfate_dihydrate SDN:P01::CCOMD002 SDN:P06::UMMC CDOM_fluor_WETLabs nan Concentration of coloured dissolved organic matter {CDOM Gelbstoff} per unit volume of the water body [dissolved plus reactive particulate phase] by in-situ WET Labs FDOM ECO fluorometer Milligrams per cubic metre mg/m^3 mg/m**3 nan
Cdom CDOM fluorescence mg/m^3 concentration_of_colored_dissolved_organic_matter_in_sea_water_expressed_as_equivalent_mass_fraction_of_quinine_sulfate_dihydrate nan SDN:P06::UMMC CDOM_fluor_WETLabs CDOM nan Milligrams per cubic metre nan nan nan
NO3 Nitrate (NO3-N) mmol/m^3 mole_concentration_of_nitrate_in_sea_water SDN:P01::NTRAZZXX SDN:P06::MMCM nan NTRA Concentration of nitrate {NO3- CAS 14797-55-8} per unit volume of the water body [unknown phase] Millimoles per cubic metre mmol/m^3 mmol/m**3 nan
Ntra Nitrate (NO3-N) mmol/m^3 mole_concentration_of_nitrate_in_sea_water SDN:P01::NTRAZZXX SDN:P06::MMCM NO3 NTRA Concentration of nitrate {NO3- CAS 14797-55-8} per unit volume of the water body [unknown phase] Millimoles per cubic metre nan nan nan
Nitrate Nitrate (NO3-N) mmol/m^3 mole_concentration_of_nitrate_in_sea_water SDN:P01::NTRAZZXX SDN:P06::MMCM NO3 NTRA Concentration of nitrate {NO3- CAS 14797-55-8} per unit volume of the water body [unknown phase] Millimoles per cubic metre nan nan nan
NTRA Nitrate (NO3-N) mmol/m^3 mole_concentration_of_nitrate_in_sea_water SDN:P01::NTRAZZXX SDN:P06::MMCM NO3 nan Concentration of nitrate {NO3- CAS 14797-55-8} per unit volume of the water body [unknown phase] Millimoles per cubic metre nan nan nan
Par Sea Water PAR uEinsteins/s/m^2 downwelling_photosynthetic_photon_flux_in_sea_water SDN:P01::IRRDUV01 SDN:P06::UMES SubsurVPAR PSAR Downwelling vector irradiance as photons of electromagnetic radiation (PAR wavelengths) in the water body by cosine-collector radiometer MicroEinsteins per square metre per second uEinsteins/s/m\^2 ueinsteins/s/m\^2 ueinsteins/s/m**2
PSAR Sea Water PAR uEinsteins/s/m^2 downwelling_photosynthetic_photon_flux_in_sea_water SDN:P01::IRRDUV01 SDN:P06::UMES SubsurVPAR nan Downwelling vector irradiance as photons of electromagnetic radiation (PAR wavelengths) in the water body by cosine-collector radiometer MicroEinsteins per square metre per second uEinsteins/s/m\^2 ueinsteins/s/m\^2 ueinsteins/s/m**2
SPar Surface PAR uEinsteins/s/m^2 surface_downwelling_photosynthetic_photon_flux_in_sea_water SDN:P01::IRRDSV01 SDN:P06::UMES SurfVPAR SPAR Downwelling vector irradiance as photons of electromagnetic radiation (PAR wavelengths) in the atmosphere by cosine-collector radiometer MicroEinsteins per square metre per second uEinsteins/s/m\^2 ueinsteins/s/m\^2 ueinsteins/s/m**2
SPAR Surface PAR uEinsteins/s/m^2 surface_downwelling_photosynthetic_photon_flux_in_sea_water SDN:P01::IRRDSV01 SDN:P06::UMES SurfVPAR nan Downwelling vector irradiance as photons of electromagnetic radiation (PAR wavelengths) in the atmosphere by cosine-collector radiometer MicroEinsteins per square metre per second uEinsteins/s/m\^2 ueinsteins/s/m\^2 ueinsteins/s/m**2
pH Sea water pH 1 sea_water_ph_reported_on_total_scale SDN:P01::PHMASSXX SDN:P06::UUPH nan nan pH (total scale) {pH[T]} of the water body pH units nan nan nan
Bottle nan nan nan nan nan BOPO nan nan nan nan nan None
PrDM Sea Pressure (Sea Surface - 0) decibars sea_water_pressure_due_to_sea_water SDN:P01::PRESPR01 SDN:P06::UPDB Pres_Z PRES Pressure (spatial coordinate) exerted by the water body by profiling pressure sensor and correction to read zero at sea level Decibars nan nan nan
T090C Temperature (ITS-90) degC sea_water_temperature SDN:P01::TEMPS901 SDN:P06::UPAA CTDTmp90 TE90 Temperature (ITS-90) of the water body by CTD or STD Degrees Celsius degC deg C degrees C
T190C Temperature (ITS-90) degC sea_water_temperature SDN:P01::TEMPS902 SDN:P06::UPAA CTDTmp90_2 TE90 Temperature (ITS-90) of the water body by CTD or STD (second sensor) Degrees Celsius degC deg C degrees C
WetCDOM CDOM fluorescence mg/m^3 concentration_of_colored_dissolved_organic_matter_in_sea_water_expressed_as_equivalent_mass_fraction_of_quinine_sulfate_dihydrate SDN:P01::CCOMD002 SDN:P06::UMMC CDOM_fluor_WETLabs CDOM Concentration of coloured dissolved organic matter {CDOM Gelbstoff} per unit volume of the water body [dissolved plus reactive particulate phase] by in-situ WET Labs FDOM ECO fluorometer Milligrams per cubic metre nan nan nan
Sal00 Practical Salinity PSU sea_water_practical_salinity SDN:P01::PSALST01 SDN:P06::PSUX P_sal_CTD PSAL Practical salinity of the water body by CTD and computation using UNESCO 1983 algorithm Practical salinity units PSU psu 1e-3
FlSP chl-a fluorescence ug/L mass_concentration_of_chlorophyll_a_in_sea_water SDN:P01::CPHLPR01 SDN:P06::UGPL chla_water_ISfluor FLOR Concentration of chlorophyll-a {chl-a CAS 479-61-8} per unit volume of the water body [particulate >unknown phase] by in-situ chlorophyll fluorometer Micrograms per litre ug/L mg/m\^3 nan
Sigma-t00 Sigma-T Rho(S, T, 0)-1000 kg/m^3 sea_water_sigma_t SDN:P01::SIGTEQST SDN:P06::UKMC SigmaT SIGT Sigma-T of the water body by computation from salinity and temperature using UNESCO algorithm Kilograms per cubic metre nan nan nan
Sigma-é00 Sigma-Theta (Rho(S, Theta, 0)-1000) kg/m^3 sea_water_sigma_theta SDN:P01::SIGTEQ01 SDN:P06::UKMC SigmaTheta SIGP Sigma-theta of the water body by computation from salinity and potential temperature using UNESCO algorithm Kilograms per cubic metre nan nan nan
Spar Surface Photosynthetic Active Radiation uEinsteins/s/m^2 surface_downwelling_photosynthetic_photon_flux_in_sea_water SDN:P01::IRRDSV01 SDN:P06::UMES SurfVPAR SPAR Downwelling vector irradiance as photons of electromagnetic radiation (PAR wavelengths) in the atmosphere by cosine-collector radiometer MicroEinsteins per square metre per second uEinsteins/s/m\^2 µeinsteins/s/m\^2 nan
Sbeox0Mm/L Dissolved Oxygen uM mole_concentration_of_dissolved_molecular_oxygen_in_sea_water SDN:P01::DOXYZZ01 SDN:P06::UPOX WC_dissO2_IS OXYM Concentration of oxygen {O2 CAS 7782-44-7} per unit volume of the water body [dissolved plus reactive particulate phase] by in-situ sensor Micromoles per litre nan nan nan
CStarAt0 Light Transmission % nan SDN:P01::POPTDR01 SDN:P06::UPCT Trans_Red_25cm nan Transmittance (red light wavelength) per 25cm of the water body by 25cm path length red light transmissometer Percent nan nan nan
Latitude Latitude nan latitude SDN:P01::ALATZZ01 SDN:P06::DEGN nan latitude Latitude north Degrees north nan nan nan
Longitude Longitude nan longitude SDN:P01::ALONZZ01 SDN:P06::DEGE Lon longitude Longitude east Degrees east nan nan nan
Upoly0 Nitrate (NO3-N) mmol/m^3 mole_concentration_of_nitrate_in_sea_water SDN:P01::NTRAZZXX SDN:P06::MMCM nan NTRA Concentration of nitrate {NO3- CAS 14797-55-8} per unit volume of the water body [unknown phase] Millimoles per cubic metre nan nan nan