-
Notifications
You must be signed in to change notification settings - Fork 13
Fix ssp correction #69
New issue
Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.
By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.
Already on GitHub? Sign in to your account
Open
krisvanneste
wants to merge
11
commits into
SeismicSource:main
Choose a base branch
from
krisvanneste:fix_ssp_correction
base: main
Could not load branches
Branch not found: {{ refName }}
Loading
Could not load tags
Nothing to show
Loading
Are you sure you want to change the base?
Some commits from the old base branch may be removed from the timeline,
and old review comments may become outdated.
Open
Changes from all commits
Commits
Show all changes
11 commits
Select commit
Hold shift + click to select a range
39ecbc0
Limit both the original spectrum and the residual to their common fre…
krisvanneste 453e1b6
Remove possible NaN values from log-spaced spectrum and associated we…
krisvanneste 2fa769c
Added modification in residual correction to changelog.
krisvanneste eb242d5
ssp_correction: move interpolation to a separate function
claudiodsf 8a539f5
Remove unnexessary return from station_correction()
claudiodsf 34a0433
Exit with error if building the weight from the noise fails
claudiodsf cbf3021
Reimplemented station_correction function without altering the freque…
krisvanneste 9d3be4d
Replace NaN with small value in weight data
claudiodsf 6053c27
Allow for NaN values in stacked spectra
claudiodsf ade7fed
Ignore NaN data in computing radiated energy
claudiodsf e818baa
Move interpolation into Spectrum class
claudiodsf File filter
Filter by extension
Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
There are no files selected for viewing
This file contains hidden or bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
|
@@ -24,6 +24,8 @@ Copyright (c) 2011-2025 Claudio Satriano <[email protected]> | |
- New config parameter `clipping_min_amplitude_ratio` to set a threshold for | ||
trace amplitude below which the trace is not checked for clipping | ||
- Use spectral interpolation to compute and apply station residuals | ||
- Limit spectrum and residual to common frequency range when applying | ||
correction before fitting | ||
|
||
### Post-Inversion | ||
|
||
|
This file contains hidden or bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
This file contains hidden or bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
|
@@ -9,6 +9,7 @@ | |
|
||
:copyright: | ||
2012-2025 Claudio Satriano <[email protected]> | ||
Kris Vanneste <[email protected]> | ||
:license: | ||
CeCILL Free Software License Agreement v2.1 | ||
(http://www.cecill.info/licences.en.html) | ||
|
@@ -22,6 +23,7 @@ | |
import math | ||
import yaml | ||
import numpy as np | ||
from scipy.interpolate import interp1d | ||
# Reduce logging loevel for h5py. | ||
# For h5py, this has to be done before importing the module. | ||
logging.getLogger('h5py').setLevel(logging.WARNING) | ||
|
@@ -83,6 +85,54 @@ def __deepcopy__(self, memo): | |
return new_dict | ||
|
||
|
||
def _find_nan_edges(data): | ||
""" | ||
Helper function to find NaN values at the beginning and end of a 1D array. | ||
|
||
:param data: The data to check for NaN values. | ||
:return: The indices of the NaN values at the beginning and end | ||
of the data. | ||
If there are no NaN values, an empty array is returned. | ||
If all values are NaN, an array with all indices is returned. | ||
""" | ||
nan_mask = np.isnan(data) | ||
# All NaNs | ||
if not np.any(~nan_mask): | ||
return np.arange(len(data)) | ||
# NaNs at the beginning | ||
start_nan_idxs = np.where(nan_mask[:np.argmax(~nan_mask)])[0] | ||
# NaNs at the end | ||
end_nan_mask_reversed = nan_mask[::-1] | ||
end_nan_idxs = np.where( | ||
end_nan_mask_reversed[:np.argmax(~end_nan_mask_reversed)])[0] | ||
end_nan_idxs = len(data) - 1 - end_nan_idxs | ||
# Concatenate the two arrays. The result can be empty if there are | ||
# no NaNs at the beginning or end. | ||
return np.concatenate((start_nan_idxs, end_nan_idxs)) | ||
|
||
|
||
def _interpolate_data_to_new_freq( | ||
data, freq, new_freq, fill_value='extrapolate', fix_negative=True): | ||
""" | ||
Helper function to interpolate data to a new frequency range. | ||
|
||
:param data: The data to interpolate. | ||
:param freq: The original frequency range. | ||
:param new_freq: The new frequency range. | ||
:param fill_value: The value to use for extrapolation. | ||
If 'extrapolate', the data is extrapolated to the new frequency range. | ||
See scipy.interpolate.interp1d for more details. | ||
:param fix_negative: If True, negative values are replaced by the | ||
minimum value of the original data. | ||
:return: The interpolated data. | ||
""" | ||
f = interp1d(freq, data, fill_value=fill_value, bounds_error=False) | ||
new_data = f(new_freq) | ||
if fix_negative: | ||
new_data[new_data <= 0] = np.min(data) | ||
return new_data | ||
|
||
|
||
class Spectrum(): | ||
""" | ||
A class to handle amplitude spectra. | ||
|
@@ -320,6 +370,176 @@ def freq_logspaced(self, value): | |
self.stats.delta_logspaced = 1. | ||
self._freq_logspaced = value | ||
|
||
def make_freq_logspaced(self, delta_logspaced=None): | ||
""" | ||
Create a logspaced frequency axis from the linear frequency axis. | ||
|
||
If logspaced data already exist, it is reinterpolated to the new | ||
logspaced frequencies. | ||
|
||
:param delta_logspaced: The log10([Hz]) sample interval. | ||
If None, it is computed from the linear frequency axis. | ||
""" | ||
if self.freq.size == 0: | ||
raise ValueError('Frequency axis is empty') | ||
if np.any(self.freq <= 0): | ||
raise ValueError('Frequency axis must be positive') | ||
log_freq = np.log10(self.freq) | ||
if delta_logspaced is None: | ||
log_df = log_freq[-1] - log_freq[-2] | ||
else: | ||
log_df = delta_logspaced | ||
freq_logspaced =\ | ||
10**(np.arange(log_freq[0], log_freq[-1] + log_df, log_df)) | ||
# If logspaced frequencies already exist, | ||
# reinterpolate the data to the new logspaced frequencies | ||
if ( | ||
self.data_logspaced.size > 0 | ||
and self.freq_logspaced.size != freq_logspaced.size | ||
): | ||
self.interp_data_logspaced_to_new_freq(freq_logspaced) | ||
self.freq_logspaced = freq_logspaced | ||
|
||
def make_logspaced_from_linear(self, which='data'): | ||
""" | ||
Convert the linear data to logspaced data. | ||
|
||
:param which: The data to convert. | ||
One of 'data', 'data_mag' or 'both'. | ||
|
||
:note: The logspaced frequency axis must exist. | ||
""" | ||
if which == 'data': | ||
data = self.data | ||
fix_negative = True | ||
elif which == 'data_mag': | ||
data = self.data_mag | ||
fix_negative = False | ||
elif which == 'both': | ||
self.make_logspaced_from_linear('data') | ||
self.make_logspaced_from_linear('data_mag') | ||
return | ||
else: | ||
raise ValueError( | ||
f'Invalid value for "which": {which}. ' | ||
'Must be one of "data", "data_mag" or "both".' | ||
) | ||
if self.freq.size == 0: | ||
raise ValueError('Frequency axis is empty') | ||
if data.size == 0: | ||
raise ValueError('Data axis is empty') | ||
if self.freq_logspaced.size == 0: | ||
raise ValueError('Logspaced frequency axis is empty') | ||
# Find NaN indexes at the beginning and end of the data | ||
nan_idxs = _find_nan_edges(data) | ||
# Do not extrapolate if there are NaN values at the beginning or end | ||
fill_value = np.nan if nan_idxs.size else 'extrapolate' | ||
data = np.delete(data, nan_idxs) | ||
freq = np.delete(self.freq, nan_idxs) | ||
# Reinterpolate data using log10 frequencies | ||
data_logspaced = _interpolate_data_to_new_freq( | ||
data, freq, self.freq_logspaced, | ||
fill_value=fill_value, fix_negative=fix_negative | ||
) | ||
if which == 'data': | ||
self.data_logspaced = data_logspaced | ||
else: | ||
self.data_mag_logspaced = data_logspaced | ||
|
||
def make_linear_from_logspaced(self, which='data_logspaced'): | ||
""" | ||
Convert the logspaced data to linear data. | ||
|
||
:param which: The data to convert. | ||
One of 'data_logspaced', 'data_mag_logspaced' or 'both'. | ||
|
||
:note: The linear frequency axis must exist. | ||
""" | ||
if which == 'data_logspaced': | ||
data_logspaced = self.data_logspaced | ||
fix_negative = True | ||
elif which == 'data_mag_logspaced': | ||
data_logspaced = self.data_mag_logspaced | ||
fix_negative = False | ||
elif which == 'both': | ||
self.make_linear_from_logspaced('data_logspaced') | ||
self.make_linear_from_logspaced('data_mag_logspaced') | ||
return | ||
else: | ||
raise ValueError( | ||
f'Invalid value for "which": {which}. ' | ||
'Must be one of "data_logspaced", "data_mag_logspaced" or ' | ||
'"both".' | ||
) | ||
if self.freq_logspaced.size == 0: | ||
raise ValueError('Logspaced frequency axis is empty') | ||
if self.data_logspaced.size == 0: | ||
raise ValueError('Data logspaced axis is empty') | ||
if self.freq.size == 0: | ||
raise ValueError('Frequency axis is empty') | ||
# Reinterpolate data using linear frequencies | ||
data = _interpolate_data_to_new_freq( | ||
data_logspaced, self.freq_logspaced, self.freq, | ||
fix_negative=fix_negative | ||
) | ||
if which == 'data_logspaced': | ||
self.data = data | ||
else: | ||
self.data_mag = data | ||
|
||
def interp_data_to_new_freq(self, new_freq, fill_value='extrapolate'): | ||
""" | ||
Interpolate the linear data to a new frequency axis. | ||
|
||
:param new_freq: The new frequency axis. | ||
:param fill_value: The value to use for extrapolation. | ||
If 'extrapolate', the data is extrapolated to the new frequency | ||
range. See scipy.interpolate.interp1d for more details. | ||
""" | ||
if self.freq.size == 0: | ||
raise ValueError('Frequency axis is empty') | ||
if self.data.size == 0: | ||
raise ValueError('Data axis is empty') | ||
self.data = _interpolate_data_to_new_freq( | ||
self.data, self.freq, new_freq, | ||
fill_value=fill_value, | ||
fix_negative=True | ||
) | ||
if self.data_mag.size: | ||
self.data_mag = _interpolate_data_to_new_freq( | ||
self.data_mag, self.freq, new_freq, | ||
fill_value=fill_value, | ||
fix_negative=False | ||
) | ||
self.freq = new_freq | ||
|
||
def interp_data_logspaced_to_new_freq( | ||
self, new_freq, fill_value='extrapolate'): | ||
""" | ||
Interpolate the logspaced data to a new frequency axis. | ||
|
||
:param new_freq: The new frequency axis. | ||
:param fill_value: The value to use for extrapolation. | ||
If 'extrapolate', the data is extrapolated to the new frequency | ||
range. See scipy.interpolate.interp1d for more details. | ||
""" | ||
if self.freq_logspaced.size == 0: | ||
raise ValueError('Logspaced frequency axis is empty') | ||
if self.data_logspaced.size == 0: | ||
raise ValueError('Data logspaced axis is empty') | ||
self.data_logspaced = _interpolate_data_to_new_freq( | ||
self.data_logspaced, self.freq_logspaced, new_freq, | ||
fill_value=fill_value, | ||
fix_negative=True | ||
) | ||
if self.data_mag_logspaced.size: | ||
self.data_mag_logspaced = _interpolate_data_to_new_freq( | ||
self.data_mag_logspaced, self.freq_logspaced, new_freq, | ||
fill_value=fill_value, | ||
fix_negative=False | ||
) | ||
self.freq_logspaced = new_freq | ||
|
||
def copy(self): | ||
"""Return a copy of the spectrum.""" | ||
spec_copy = Spectrum() | ||
|
This file contains hidden or bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Oops, something went wrong.
Add this suggestion to a batch that can be applied as a single commit.
This suggestion is invalid because no changes were made to the code.
Suggestions cannot be applied while the pull request is closed.
Suggestions cannot be applied while viewing a subset of changes.
Only one suggestion per line can be applied in a batch.
Add this suggestion to a batch that can be applied as a single commit.
Applying suggestions on deleted lines is not supported.
You must change the existing code in this line in order to create a valid suggestion.
Outdated suggestions cannot be applied.
This suggestion has been applied or marked resolved.
Suggestions cannot be applied from pending reviews.
Suggestions cannot be applied on multi-line comments.
Suggestions cannot be applied while the pull request is queued to merge.
Suggestion cannot be applied right now. Please check back later.
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
Change to
nan_idxs = []
if you want to ignore this step.There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
Claudio,
Thanks for the update. I will try to run it one of the coming days (Friday at the latest). Is that OK?
Kris
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
Sure! thanks!