→ Plot PSD-RMS
PSD-RMS Visualisation
Plots derived from the psd_rms workflow step: time-series, clock plots,
hour-maps, grid-maps, and daily stacks. Ported and adapted from the
SeismoRMS / SeismoSocialDistancing project
(https://github.com/ThomasLecocq/SeismoRMS).
CLI usage:
msnoise qc plot_psd_rms BE.UCC..HHZ --type timeseries --band 1.0-10.0
msnoise qc plot_psd_rms BE.UCC..HHZ --type clockplot --timezone Europe/Brussels
msnoise qc plot_psd_rms BE.UCC..HHZ --type hourmap
msnoise qc plot_psd_rms BE.UCC..HHZ --type gridmap
msnoise qc plot_psd_rms BE.UCC..HHZ --type dailyplot
Notebook usage:
from msnoise.core.db import connect
from msnoise.results import MSNoiseResult
from msnoise.plots.psd_rms import load_rms, plot_timeseries, plot_clockplot
db = connect()
r = MSNoiseResult.from_ids(db, psd=1, psd_rms=1)
data = load_rms(r, ["BE.UCC..HHZ"])
plot_timeseries(data, band="1.0-10.0")
- class msnoise.plots.psd_rms.DayNightWindow(day_start: float = 7.0, day_end: float = 19.0, day_color: str = '#fffacd', night_color: str = '#cce5ff', day_alpha: float = 0.15, night_alpha: float = 0.12)
Define the local “daytime” and “nighttime” hour windows.
Used by all plot functions to:
draw the daytime median overlay on time-series plots,
shade night hours on clockplots, gridmaps, and dailyplots,
shade night sectors on hourmaps.
Hours are expressed as floats in 0..24 in local time (after timezone conversion). The window wraps midnight correctly, so day_start=22 with day_end=6 defines a “daytime” from 22:00 to 06:00 (e.g. a city station where the quiet hours are at night).
- Parameters:
- day_start:
Start of daytime in decimal hours (default
7.0= 07:00 local).- day_end:
End of daytime in decimal hours (default
19.0= 19:00 local).- day_color:
Matplotlib colour for daytime shading (default light yellow).
- night_color:
Matplotlib colour for nighttime shading (default light blue).
- day_alpha:
Alpha for daytime shading (default
0.15).- night_alpha:
Alpha for nighttime shading (default
0.12).
Examples
>>> # Standard business-hours window >>> dnw = DayNightWindow(day_start=8, day_end=18)
>>> # Seismological quiet window: anthropogenic noise lowest 22:00-05:00 >>> dnw = DayNightWindow(day_start=5, day_end=22)
>>> # Marine station: no meaningful day/night -- disable shading >>> dnw = DayNightWindow(day_start=0, day_end=24)
- shade_cartesian_hour(ax) None
Shade day / night horizontal bands on a Cartesian hour-of-day y-axis.
Assumes the y-axis represents hours 0-24 (gridmap).
- msnoise.plots.psd_rms.load_rms(result, seed_ids: list[str] | None = None, *, psd_id: int = 1, psd_rms_id: int = 1) dict[str, pandas.core.frame.DataFrame]
Load PSD-RMS results via an
MSNoiseResult.This is the preferred entry point for notebooks and scripts. It delegates all path resolution to
MSNoiseResult.get_psd_rms(), which correctly handles thepsd_N/psd_rms_N/_output/lineage without hard-coding step names or output folders.- Parameters:
- result:
An
MSNoiseResultthat containspsd_rmsin its lineage, obtained via:from msnoise.results import MSNoiseResult r = MSNoiseResult.from_ids(db, psd=1, psd_rms=1)
- seed_ids:
Optional list of SEED IDs to load (e.g.
["BE.UCC..HHZ"]). WhenNone(default) every station found on disk is returned.
- Returns:
- dict
Mapping seed_id ->
DataFrame(columns = band labels, index =DatetimeIndexin UTC). Stations with no data on disk are silently omitted.
Examples
Preferred – MSNoiseResult:
>>> from msnoise.core.db import connect >>> from msnoise.results import MSNoiseResult >>> from msnoise.plots.psd_rms import load_rms >>> db = connect() >>> r = MSNoiseResult.from_ids(db, psd=1, psd_rms=1) >>> data = load_rms(r) # all stations auto-discovered >>> data = load_rms(r, ["BE.UCC..HHZ"]) # single station
- msnoise.plots.psd_rms.plot_timeseries(data: dict[str, pandas.core.frame.DataFrame], band: str | None = None, scale: float = 1000000000.0, unit: str = 'nm', annotations: dict | None = None, time_zone: str = 'UTC', resample_freq: str = '30min', agg_func: str | Callable = 'mean', day_night: DayNightWindow | None = None, logo: str | None = None, show: bool = True, outfile: str | None = None) Figure
Time-series plot of RMS amplitudes for one or multiple stations.
- Parameters:
- data:
dict
{seed_id: DataFrame}as returned byload_rms().- band:
Frequency band label (column name in the DataFrame, e.g.
"1.0-10.0"). Defaults to the first available band.- scale:
Multiply amplitudes before display (default
1e9-> nm for DISP).- unit:
Y-axis unit label (default
"nm").- annotations:
Optional dict
{date_string: label}for vertical event lines.- time_zone:
Time zone for x-axis localisation (default
"UTC").- resample_freq:
Resampling frequency string (default
"30min").- agg_func:
Aggregation function applied to each resampled bin. Accepts a string (
"mean","median","max", …) or any callable (e.g.np.nanmedian,lambda s: s.quantile(0.1)). Default"mean".- day_night:
DayNightWindowdefining local daytime/nighttime hours. Controls: (1) the daytime median overlay time window and label, and (2) night shading applied as vertical spans. Defaults toDEFAULT_DAY_NIGHT(07:00-19:00). PassDayNightWindow(day_start=0, day_end=24)to disable shading.- logo:
URL or path to a logo image to embed bottom-left (optional).
- show:
Call
plt.show()when done.- outfile:
Save figure to this path.
- Returns:
- msnoise.plots.psd_rms.plot_clockplot(data: dict[str, pandas.core.frame.DataFrame], band: str | None = None, scale: float = 1000000000.0, unit: str = 'nm', split_date: str | None = None, time_zone: str = 'UTC', resample_freq: str = '30min', agg_func: str | Callable = 'mean', day_night: DayNightWindow | None = None, show: bool = True, outfile: str | None = None) Figure
Polar 24-hour clock plot of median RMS per weekday / hour.
When split_date is provided two panels are drawn (before / after), useful for comparing noise regimes across an event or intervention.
- Parameters:
- data:
dict
{seed_id: DataFrame}as returned byload_rms(). Only the first station is used; for multi-station comparisons build a side-by-side figure usingstack_wday_time()directly (see the multi-station notebook example).- band:
Frequency band label.
- scale:
Amplitude scale factor (default
1e9-> nm).- unit:
Unit label for radial tick labels.
- split_date:
ISO date string (
"YYYY-MM-DD") to split the data. WhenNonea single panel is drawn covering the full time range.- time_zone:
Local time zone for hour-of-day grouping (default
"UTC").- resample_freq:
Resampling frequency string (default
"30min").- agg_func:
Aggregation function for resampled bins (
"mean","median", callable, …). Default"mean".- day_night:
DayNightWindowfor shading night sectors on the polar axes. Defaults toDEFAULT_DAY_NIGHT(07:00-19:00).- show:
Call
plt.show().- outfile:
Save figure to this path.
- Returns:
- msnoise.plots.psd_rms.plot_hourmap(data: dict[str, pandas.core.frame.DataFrame], band: str | None = None, scale: float = 1000000000.0, unit: str = 'nm', annotations: dict | None = None, time_zone: str = 'UTC', resample_freq: str = '30min', agg_func: str | Callable = 'mean', day_night: DayNightWindow | None = None, show: bool = True, outfile: str | None = None) Figure
Polar hour-map: radial axis = elapsed days, angular axis = hour of day.
Equivalent to SeismoRMS
clockmap/hourmap. Best for datasets up to ~1-2 years; beyond that the polar readability degrades andplot_gridmap()is preferable.- Parameters:
- data:
dict
{seed_id: DataFrame}as returned byload_rms().- band:
Frequency band label.
- scale / unit:
Amplitude scale and unit label.
- annotations:
Optional dict
{date_str: label}to mark events as radial markers.- time_zone:
Local time zone.
- resample_freq:
Resampling frequency.
- agg_func:
Aggregation function for resampled bins (
"mean","median", callable, …). Default"mean".- day_night:
DayNightWindowfor shading night sectors on the polar axes. Defaults toDEFAULT_DAY_NIGHT(07:00-19:00).- show / outfile:
Display / save options.
- Returns:
- msnoise.plots.psd_rms.plot_gridmap(data: dict[str, pandas.core.frame.DataFrame], band: str | None = None, scale: float = 1000000000.0, unit: str = 'nm', annotations: dict | None = None, time_zone: str = 'UTC', resample_freq: str = '30min', agg_func: str | Callable = 'mean', day_night: DayNightWindow | None = None, show: bool = True, outfile: str | None = None) Figure
Cartesian grid-map: x = calendar date, y = hour of day.
Equivalent to SeismoRMS
gridmap. Handles long time series (years) better thanplot_hourmap().- Parameters:
- data:
dict
{seed_id: DataFrame}as returned byload_rms().- band:
Frequency band label.
- scale / unit:
Amplitude scale and unit label.
- annotations:
Optional dict
{date_str: label}for event markers.- time_zone:
Local time zone.
- resample_freq:
Resampling frequency.
- agg_func:
Aggregation function for resampled bins (
"mean","median", callable, …). Default"mean".- day_night:
DayNightWindowfor shading night bands on the hour-of-day y-axis. Defaults toDEFAULT_DAY_NIGHT(07:00-19:00).- show / outfile:
Display / save options.
- Returns:
- msnoise.plots.psd_rms.plot_dailyplot(data: dict[str, pandas.core.frame.DataFrame], band: str | None = None, scale: float = 1000000000.0, unit: str = 'nm', split_date: str | None = None, time_zone: str = 'UTC', resample_freq: str = '30min', agg_func: str | Callable = 'mean', day_night: DayNightWindow | None = None, show: bool = True, outfile: str | None = None) Figure
Hour-of-day stacked plot, one curve per weekday.
Solid curves cover the full dataset (or the period before split_date when provided); dashed curves represent the period after split_date.
- Parameters:
- data:
dict
{seed_id: DataFrame}as returned byload_rms().- band:
Frequency band label.
- scale / unit:
Amplitude scale and unit label.
- split_date:
Optional date string to draw a pre/post comparison (solid vs dashed).
- time_zone:
Local time zone for hour-of-day grouping.
- resample_freq:
Resampling frequency.
- agg_func:
Aggregation function for resampled bins (
"mean","median", callable, …). Default"mean".- day_night:
DayNightWindowfor shading night bands on the hour-of-day x-axis. Defaults toDEFAULT_DAY_NIGHT(07:00-19:00).- show / outfile:
Display / save options.
- Returns:
- msnoise.plots.psd_rms.main(seed_ids: list[str], psd_id: int = 1, psd_rms_id: int = 1, plot_type: str = 'timeseries', band: str | None = None, scale: float = 1000000000.0, unit: str = 'nm', time_zone: str = 'UTC', day_start: float = 7.0, day_end: float = 19.0, split_date: str | None = None, annotations: dict | None = None, outfile: str | None = None, show: bool = True, loglevel: str = 'INFO') None
CLI/script entry point: load RMS data and dispatch to plot function(s).
Constructs an
MSNoiseResultfrom the database, loads data viaload_rms(), then dispatches to the requested plot function(s).- Parameters:
- seed_ids:
List of SEED IDs (e.g.
["BE.UCC..HHZ"]). Empty list -> all stations found on disk for this lineage.- psd_id:
Config-set number for the
psdstep (default1).- psd_rms_id:
Config-set number for the
psd_rmsstep (default1).- plot_type:
One of
"timeseries","clockplot","hourmap","gridmap","dailyplot", or"all".- band:
Frequency band column label (e.g.
"1.0-10.0").- scale / unit:
Amplitude scale and label.
- time_zone:
Time zone string for local-time plots.
- day_start:
Start of daytime in decimal hours (default
7.0).- day_end:
End of daytime in decimal hours (default
19.0).- split_date:
Date string for before/after split (clockplot, dailyplot).
- annotations:
Dict
{date_str: label}for event markers (timeseries, hourmap, gridmap).- outfile:
Base filename for saved figures; a
_<type>suffix is inserted before the extension for each type whenplot_type="all".- show:
Call
plt.show().- loglevel:
Logging level string.
- msnoise.plots.psd_rms.localize_tz_and_reindex(df: DataFrame, freq: str = '15min', time_zone: str = 'UTC', agg_func: str | Callable = 'mean') DataFrame
Convert a UTC-indexed DataFrame to time_zone, resample to freq, and aggregate with agg_func.
- Parameters:
- df:
DataFrame with a timezone-naive UTC
DatetimeIndex.- freq:
Pandas offset string for resampling (default
"15min").- time_zone:
Target timezone string (default
"UTC").- agg_func:
Aggregation function applied after resampling. Can be:
a string supported by
agg()("mean","median","max","min","std", …)any callable that accepts a
Seriesand returns a scalar (e.g.np.nanmedian,lambda s: s.quantile(0.1)).
Default
"mean"preserves the original behaviour.
- Returns:
- msnoise.plots.psd_rms.pivot_for_hourmap(data: DataFrame, columns: str = 'angles') DataFrame
Pivot a single-column RMS DataFrame into a (days x time) matrix.
- Parameters:
- data:
DataFrame with a
DatetimeIndexand a single RMS column (already scaled to the desired display unit).- columns:
"angles"(default) converts the time axis to radians 0-2pi for polar pcolormesh;"hours"leaves them as float hours 0-24 for the Cartesian grid-map.
- Returns:
DataFrameRows = elapsed integer days since the first sample; columns = hour or angle values.
- msnoise.plots.psd_rms.stack_wday_time(df: DataFrame, scale: float = 1.0) DataFrame
Median-stack a single-column RMS DataFrame into a (hour x weekday) table.
- Parameters:
- df:
DataFrame with a
DatetimeIndexand a single column. The column name is ignored — only the first column is used.- scale:
Multiply all values by this factor (default
1.0).
- Returns:
DataFrameIndex = float hour-of-day; columns = weekday names ordered Monday to Sunday.