Skip to content
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
1 change: 1 addition & 0 deletions changelog/158.feature.rst
Original file line number Diff line number Diff line change
@@ -0,0 +1 @@
Restored the WIND/WAVES ``Fido`` client `~radiospectra.net.sources.wind.WAVESClient` and updated it to use NASA SPDF server and altered paths and filenames.
2 changes: 2 additions & 0 deletions radiospectra/net/__init__.py
Original file line number Diff line number Diff line change
Expand Up @@ -4,11 +4,13 @@
from radiospectra.net.sources.ilofar import ILOFARMode357Client
from radiospectra.net.sources.psp import RFSClient
from radiospectra.net.sources.rstn import RSTNClient
from radiospectra.net.sources.wind import WAVESClient

__all__ = [
"eCALLISTOClient",
"EOVSAClient",
"RFSClient",
"RSTNClient",
"WAVESClient",
"ILOFARMode357Client",
]
99 changes: 99 additions & 0 deletions radiospectra/net/sources/tests/test_wind_client.py
Original file line number Diff line number Diff line change
@@ -0,0 +1,99 @@
from datetime import datetime
from unittest import mock

import numpy as np
import pytest

import astropy.units as u

from sunpy.net import Fido
from sunpy.net import attrs as a

from radiospectra.net.sources.wind import WAVESClient

MOCK_PATH = "sunpy.net.scraper.urlopen"


@pytest.fixture
def client():
return WAVESClient()


@pytest.mark.remote_data
def test_fido():
atr = a.Time("2020/01/02", "2020/01/03")
res = Fido.search(atr, a.Instrument("WAVES"))

assert isinstance(res[0].client, WAVESClient)
assert len(res[0]) == 4


# Taken from https://spdf.gsfc.nasa.gov/pub/data/wind/waves/rad1_idl_binary/2020/
http_cont_rad1 = """
<a href="wind_waves_rad1_20200102.R1">wind_waves_rad1_20200102.R1</a>
<a href="wind_waves_rad1_20200103.R1">wind_waves_rad1_20200103.R1</a>
"""

# Taken from https://spdf.gsfc.nasa.gov/pub/data/wind/waves/rad2_idl_binary/2020/
http_cont_rad2 = """
<a href="wind_waves_rad2_20200102.R2">wind_waves_rad2_20200102.R2</a>
<a href="wind_waves_rad2_20200103.R2">wind_waves_rad2_20200103.R2</a>
"""


@pytest.fixture
def html_responses():
return [http_cont_rad1, http_cont_rad2]


@mock.patch(MOCK_PATH)
def test_waves_client(mock_urlopen, client, html_responses):
mock_urlopen.return_value.read = mock.MagicMock()
mock_urlopen.return_value.read.side_effect = html_responses
mock_urlopen.close = mock.MagicMock(return_value=None)
atr = a.Time("2020/01/02", "2020/01/03")
query = client.search(atr)

called_urls = [
"https://spdf.gsfc.nasa.gov/pub/data/wind/waves/rad1_idl_binary/2020/",
"https://spdf.gsfc.nasa.gov/pub/data/wind/waves/rad2_idl_binary/2020/",
]
assert called_urls == [call[0][0] for call in mock_urlopen.call_args_list]
assert len(query) == 4
assert query[0]["Source"] == "WIND"
assert query[0]["Provider"] == "SPDF"

wave = [20, 1040] * u.kHz
assert np.array_equal(query[0]["Wavelength"], wave)
assert query[0]["Start Time"].datetime == datetime(2020, 1, 2)
assert query[0]["End Time"].datetime == datetime(2020, 1, 2, 23, 59, 59, 999000)

wave = [1075, 13825] * u.kHz
assert np.array_equal(query[3]["Wavelength"], wave)
assert query[3]["Start Time"].datetime == datetime(2020, 1, 3)
assert query[3]["End Time"].datetime == datetime(2020, 1, 3, 23, 59, 59, 999000)

query_urls = [row["url"] for row in query]
assert (
"https://spdf.gsfc.nasa.gov/pub/data/wind/waves/rad1_idl_binary/2020/wind_waves_rad1_20200102.R1" in query_urls
)
assert (
"https://spdf.gsfc.nasa.gov/pub/data/wind/waves/rad2_idl_binary/2020/wind_waves_rad2_20200103.R2" in query_urls
)


@pytest.mark.parametrize(
("query_wave", "receivers"),
[
(a.Wavelength(1 * u.GHz, 2 * u.GHz), []),
(a.Wavelength(1 * u.Hz, 2 * u.Hz), []),
(a.Wavelength(20 * u.kHz, 150 * u.kHz), ["rad1"]),
(a.Wavelength(1.5 * u.MHz, 15 * u.MHz), ["rad2"]),
(a.Wavelength(5 * u.MHz, 10 * u.MHz), ["rad2"]),
(a.Wavelength(100 * u.Hz, 100 * u.kHz), ["rad1"]),
(a.Wavelength(20 * u.kHz, 15 * u.MHz), ["rad1", "rad2"]),
(a.Wavelength(5 * u.kHz, 20 * u.MHz), ["rad1", "rad2"]),
],
)
def test_check_wavelength(query_wave, receivers, client):
assert set(client._check_wavelengths(query_wave)) == set(receivers)
145 changes: 145 additions & 0 deletions radiospectra/net/sources/wind.py
Original file line number Diff line number Diff line change
@@ -0,0 +1,145 @@
import astropy.units as u

from sunpy.net import attrs as a
from sunpy.net.dataretriever.client import GenericClient, QueryResponse
from sunpy.net.scraper import Scraper
from sunpy.time.timerange import TimeRange

__all__ = ["WAVESClient"]

RECEIVER_FREQUENCIES = {
"rad1": a.Wavelength(20 * u.kHz, 1040 * u.kHz),
"rad2": a.Wavelength(1.075 * u.MHz, 13.825 * u.MHz),
}

RECEIVER_EXTENSIONS = {
"rad1": "R1",
"rad2": "R2",
}


class WAVESClient(GenericClient):
"""
Provides access to WIND/WAVES IDL binary data hosted at
`NASA Goddard Space Physics Data Facility (SPDF)
<https://spdf.gsfc.nasa.gov/pub/data/wind/waves/>`__.

Examples
--------
>>> import radiospectra.net
>>> from sunpy.net import Fido, attrs as a
>>> results = Fido.search(a.Time("2020/01/01", "2020/01/02"),
... a.Instrument("WAVES")) # doctest: +REMOTE_DATA
>>> results # doctest: +REMOTE_DATA
<sunpy.net.fido_factory.UnifiedResponse object at ...>
Results from 1 Provider:
<BLANKLINE>
4 Results from the WAVESClient:
<BLANKLINE>
Start Time End Time Instrument Source Provider Wavelength
kHz
----------------------- ----------------------- ---------- ------ -------- -----------------
2020-01-01 00:00:00.000 2020-01-01 23:59:59.999 WAVES WIND SPDF 20.0 .. 1040.0
2020-01-02 00:00:00.000 2020-01-02 23:59:59.999 WAVES WIND SPDF 20.0 .. 1040.0
2020-01-01 00:00:00.000 2020-01-01 23:59:59.999 WAVES WIND SPDF 1075.0 .. 13825.0
2020-01-02 00:00:00.000 2020-01-02 23:59:59.999 WAVES WIND SPDF 1075.0 .. 13825.0
<BLANKLINE>
<BLANKLINE>
"""

pattern = (
r"https://spdf.gsfc.nasa.gov/pub/data/wind/waves/{receiver}_idl_binary/{year_path}/"
r"wind_waves_{receiver}_{{year:4d}}{{month:2d}}{{day:2d}}.{ext}"
)

@classmethod
def _check_wavelengths(cls, wavelength):
"""
Check for overlap between given wavelength and receiver frequency coverage
defined in ``RECEIVER_FREQUENCIES``.

Parameters
----------
wavelength : `sunpy.net.attrs.Wavelength`
Input wavelength range to check

Returns
-------
`list`
List of receivers names or empty list if no overlap
"""
# Input wavelength range is completely contained in one receiver range
receivers = [k for k, v in RECEIVER_FREQUENCIES.items() if wavelength in v]
# If not defined need to continue
if not receivers:
# Overlaps but not contained in, either max in low-frequency or min in high-frequency receiver
if wavelength.min in RECEIVER_FREQUENCIES["rad2"] or wavelength.max in RECEIVER_FREQUENCIES["rad2"]:
receivers.append("rad2")
if wavelength.min in RECEIVER_FREQUENCIES["rad1"] or wavelength.max in RECEIVER_FREQUENCIES["rad1"]:
receivers.append("rad1")
# min in rad1 and max in rad2
# min and max of combined rad1 and rad2 contained in given wavelength range
if a.Wavelength(RECEIVER_FREQUENCIES["rad1"].min, RECEIVER_FREQUENCIES["rad2"].max) in wavelength:
receivers = ["rad1", "rad2"]
return receivers

def search(self, *args, **kwargs):
"""
Query this client for a list of results.

Parameters
----------
*args: `tuple`
`sunpy.net.attrs` objects representing the query.
**kwargs: `dict`
Any extra keywords to refine the search.

Returns
-------
A `QueryResponse` instance containing the query result.
"""
matchdict = self._get_match_dict(*args, **kwargs)
req_wave = matchdict.get("Wavelength", None)
receivers = RECEIVER_FREQUENCIES.keys()
if req_wave is not None:
receivers = self._check_wavelengths(req_wave)

metalist = []
start_year = matchdict["Start Time"].datetime.year
end_year = matchdict["End Time"].datetime.year
tr = TimeRange(matchdict["Start Time"], matchdict["End Time"])
for receiver in receivers:
for year in range(start_year, end_year + 1):
pattern = (
self.pattern.replace("{receiver}", receiver)
.replace("{ext}", RECEIVER_EXTENSIONS[receiver])
.replace("{year_path}", str(year))
)
scraper = Scraper(format=pattern)
filesmeta = scraper._extract_files_meta(tr)
for i in filesmeta:
i["receiver"] = receiver
rowdict = self.post_search_hook(i, matchdict)
metalist.append(rowdict)

return QueryResponse(metalist, client=self)

def post_search_hook(self, exdict, matchdict):
"""
Convert receiver metadata to the receiver frequency ranges.
"""
rowdict = super().post_search_hook(exdict, matchdict)
receiver = rowdict.pop("receiver")
fr = RECEIVER_FREQUENCIES[receiver]
rowdict["Wavelength"] = u.Quantity([float(fr.min.value), float(fr.max.value)], unit=fr.unit)
return rowdict

@classmethod
def register_values(cls):
adict = {
a.Instrument: [("WAVES", "WIND - WAVES")],
a.Source: [("WIND", "WIND")],
a.Provider: [("SPDF", "NASA Goddard Space Physics Data Facility")],
a.Wavelength: [("*")],
}
return adict
11 changes: 11 additions & 0 deletions radiospectra/spectrogram/sources/tests/test_waves.py
Original file line number Diff line number Diff line change
Expand Up @@ -11,6 +11,7 @@

from radiospectra.spectrogram import Spectrogram
from radiospectra.spectrogram.sources import WAVESSpectrogram
from radiospectra.spectrogram.spectrogram_factory import SpectrogramFactory


@mock.patch("radiospectra.spectrogram.spectrogram_factory.parse_path")
Expand Down Expand Up @@ -63,3 +64,13 @@ def test_waves_rad2(parse_path_moc):
assert spec.end_time.datetime == datetime(2020, 11, 28, 23, 59)
assert spec.wavelength.min == 1.075 * u.MHz
assert spec.wavelength.max == 13.825 * u.MHz


@mock.patch("radiospectra.spectrogram.spectrogram_factory.readsav")
def test_waves_prefixed_filename_parses_date(readsav_mock):
data_array = np.zeros((256, 1441))
readsav_mock.return_value = {"arrayb": data_array}

_, meta = SpectrogramFactory._read_idl_sav(Path("wind_waves_rad1_20200711.R1"), instrument="waves")

assert meta["start_time"].isot == "2020-07-11T00:00:00.000"
2 changes: 1 addition & 1 deletion radiospectra/spectrogram/spectrogram_factory.py
Original file line number Diff line number Diff line change
Expand Up @@ -680,7 +680,7 @@ def _read_idl_sav(file, instrument=None):
# bg which is already subtracted from data ?
bg = data_array[:, -1]
data = data_array[:, :-1]
start_time = Time.strptime(file.stem, "%Y%m%d")
start_time = Time.strptime(file.stem.split("_")[-1], "%Y%m%d")
end_time = start_time + 86399 * u.s
times = start_time + (np.arange(1440) * 60 + 30) * u.s
meta = {
Expand Down