Skip to content
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

Fix SLSTR angle interpolation. #3017

Open
wants to merge 1 commit into
base: main
Choose a base branch
from
Open
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
19 changes: 16 additions & 3 deletions satpy/readers/slstr_l1b.py
Original file line number Diff line number Diff line change
Expand Up @@ -305,11 +305,24 @@
variable = variable.fillna(0)
variable.attrs["resolution"] = key.get("resolution", 1000)

from scipy.interpolate import RectBivariateSpline
spl = RectBivariateSpline(
tie_y, tie_x, variable.data[:, ::-1])

values = spl.ev(full_y, full_x)
# Check if we are interpolating angles
if "angle" in key["name"]:
# If we are interpolating the angles, we need to do so with the sine and cosine to prevent
# interpolation artifacts such as values <0 or >360.
indat = variable.data[:, ::-1]

Check notice on line 314 in satpy/readers/slstr_l1b.py

View check run for this annotation

codefactor.io / CodeFactor

satpy/readers/slstr_l1b.py#L314

Multiple spaces after operator. (E222)
sin_angles = np.sin(np.radians(indat))
cos_angles = np.cos(np.radians(indat))
sin_interp = RectBivariateSpline(tie_y, tie_x, sin_angles)
cos_interp = RectBivariateSpline(tie_y, tie_x, cos_angles)
values_sin = sin_interp.ev(full_y, full_x)
values_cos = cos_interp.ev(full_y, full_x)
values = np.degrees(np.arctan2(values_sin, values_cos)) % 360

Check warning on line 321 in satpy/readers/slstr_l1b.py

View check run for this annotation

Codecov / codecov/patch

satpy/readers/slstr_l1b.py#L314-L321

Added lines #L314 - L321 were not covered by tests
else:
# Otherwise, interpolate as normal.
spl = RectBivariateSpline(tie_y, tie_x, variable.data[:, ::-1])
values = spl.ev(full_y, full_x)

Check warning on line 325 in satpy/readers/slstr_l1b.py

View check run for this annotation

CodeScene Delta Analysis / CodeScene Cloud Delta Analysis (main)

❌ New issue: Complex Method

NCSLSTRAngles.get_dataset has a cyclomatic complexity of 10, threshold = 9. This function has many conditional statements (e.g. if, for, while), leading to lower code health. Avoid adding more conditionals and code to it without refactoring.

variable = xr.DataArray(da.from_array(values, chunks=(CHUNK_SIZE, CHUNK_SIZE)),
dims=["y", "x"], attrs=variable.attrs)
Expand Down
Loading