"""Read Datawell buoy files."""
from xarray.backends import BackendEntrypoint
from functools import cached_property
from pathlib import Path
import copy
import glob
import getpass
from datetime import datetime, timezone
from dateutil.parser import parse
import numpy as np
import xarray as xr
import pandas as pd
import wavespectra
from wavespectra.core.attributes import set_spec_attributes, attrs
from wavespectra.construct.direction import cartwright
VARIABLES = {
"tn": {"standard_name": "index", "long_name": "Transmission index"},
"tz": {
"standard_name": "sea_surface_wave_zero_upcrossing_period",
"long_name": "Zero-upcrossing period",
"units": "s",
},
"smax": {"long_name": "Maximum of the psd S(f)", "units": "m2/Hz"},
"tref": {
"standard_name": "air_temperature",
"long_name": "Reference temperature",
"units": "degC",
},
"tsea": {
"standard_name": "sea_surface_temperature",
"long_name": "Sea surface temperature",
"units": "degC",
},
"batery": {"long_name": "Battery status", "units": ""},
"av": {"long_name": "Offset of the vertical accelerometer", "units": ""},
"ax": {"long_name": "Offset of the x-accelerometer", "units": ""},
"ay": {"long_name": "Offset of the y-accelerometer", "units": ""},
"ori": {"long_name": "Buoy orientation", "units": "degree"},
"incli": {"long_name": "Magnetic inclination", "units": "degree"},
}
class Datawell:
"""Datawell SPT buoy object."""
def __init__(self, filename: str, lon=None, lat=None):
self.filename = Path(filename)
self.lon = lon
self.lat = lat
@cached_property
def location(self):
return self.filename.stem.split("}")[0]
@cached_property
def time(self):
return parse(self.filename.stem.split("}")[1]).replace(tzinfo=None)
def _set_attributes(self, dset: xr.Dataset) -> xr.Dataset:
set_spec_attributes(dset)
dset.attrs = {
"title": f"Datawell buoy spectral data from {self.filename.name}",
"source": "Datawell wave buoy",
"location": self.location,
"history": f"Generated by wavespectra v{wavespectra.__version__}",
"date_created": f"{datetime.now(timezone.utc)}",
"creator_name": getpass.getuser(),
"references": "https://www.coastalmonitoring.org/realtimedata/datawell_technicalnote_fileformats_2013-03-14.pdf",
}
# Extra non-waves attributes
for key, value in VARIABLES.items():
if key in dset:
dset[key].attrs = value
return dset
def read(self, dd: float = None) -> xr.Dataset:
dset = xr.Dataset()
with open(self.filename) as f:
# Buoy info
dset["tn"] = int(f.readline())
dset["hs"] = float(f.readline()) / 100
dset["tz"] = float(f.readline())
dset["smax"] = float(f.readline())
dset["tref"] = float(f.readline())
dset["tsea"] = float(f.readline())
dset["bat"] = int(f.readline())
dset["av"] = float(f.readline())
dset["ax"] = float(f.readline())
dset["ay"] = float(f.readline())
dset["ori"] = float(f.readline())
dset["incli"] = float(f.readline())
# Spectral data
data = pd.read_csv(f, header=None)
dset["freq"] = xr.DataArray(data[0], coords=dict(freq=data[0]), dims=("freq"))
dset["efth"] = xr.DataArray(data[1], coords=dict(freq=dset.freq)) * dset.smax
dset["dmf"] = xr.DataArray(data[2], coords=dict(freq=dset.freq))
dset["dsprf"] = xr.DataArray(data[3], coords=dict(freq=dset.freq))
dset["skewf"] = xr.DataArray(data[4], coords=dict(freq=dset.freq))
dset["kurtf"] = xr.DataArray(data[5], coords=dict(freq=dset.freq))
# 2D spectra
if dd is not None:
dir = np.arange(0, 360, dd)
dir = xr.DataArray(dir, coords=dict(dir=dir), name="dir")
cos2 = cartwright(dir=dir, dm=dset.dmf, dspr=dset.dsprf)
dset["efth"] = dset.efth * cos2
# Time dimension
dset = dset.expand_dims({"time": [self.time]}, axis=0)
# Coordinates
if self.lon is not None:
dset["lon"] = self.lon
if self.lat is not None:
dset["lat"] = self.lat
return self._set_attributes(dset)
[docs]
def read_datawell(filename_or_fileglob, dd=5.0, lon=None, lat=None) -> xr.Dataset:
"""Read Spectra from Datawell spt file.
Args:
- filename_or_fileglob (list, str): File name or file glob specifying spotter
files to read.
- dd (float): Directional spacing if 2D spectra are desired, use None to read
1D spectra.
- lon (float): Longitude of buoy.
- lat (float): Latitude of buoy.
Returns:
- dset (SpecDataset): spectra dataset object read from file.
"""
# Ensure a list of files
filenames = copy.deepcopy(filename_or_fileglob)
if isinstance(filename_or_fileglob, list):
filenames = filename_or_fileglob
elif Path(filenames).is_file():
filenames = [filename_or_fileglob]
else:
filenames = sorted(glob.glob(filename_or_fileglob))
if not filenames:
raise ValueError(f"No files found for '{filename_or_fileglob}'")
# Read files
dslist = []
for filename in filenames:
dslist.append(Datawell(filename, lon=lon, lat=lat).read(dd=dd))
dset = xr.concat(dslist, dim="time").sortby("time")
# Slice lon and lat so they are not a function of time
if lon is not None:
dset[attrs.LONNAME] = dset[attrs.LONNAME].isel(time=0, drop=True)
if lat is not None:
dset[attrs.LATNAME] = dset[attrs.LATNAME].isel(time=0, drop=True)
return dset
class DatawellBackendEntrypoint(BackendEntrypoint):
"""Datawell backend engine."""
def open_dataset(
self,
filename_or_obj,
*,
drop_variables=None,
filetype=None,
):
return read_datawell(filename_or_obj, filetype=filetype)
def guess_can_open(self, filename_or_obj):
return False
description = "Open Datawell spectra files as a wavespectra dataset."
url = "https://github.com/wavespectra/wavespectra"