eppaurora.models.ssusiq2023.ssusiq2023()   C
last analyzed

Complexity

Conditions 10

Size

Total Lines 142
Code Lines 67

Duplication

Lines 0
Ratio 0 %

Importance

Changes 0
Metric Value
cc 10
eloc 67
nop 8
dl 0
loc 142
rs 5.2798
c 0
b 0
f 0

How to fix   Long Method    Complexity    Many Parameters   

Long Method

Small methods make your code easier to understand, in particular if combined with a good name. Besides, if your method is small, finding a good name is usually much easier.

For example, if you find yourself adding comments to a method's body, this is usually a good sign to extract the commented part to a new method, and use the comment as a starting point when coming up with a good name for this new method.

Commonly applied refactorings include:

Complexity

Complex classes like eppaurora.models.ssusiq2023.ssusiq2023() often do a lot of different things. To break such a class down, we need to identify a cohesive component within that class. A common approach to find such a component is to look for fields/methods that share the same prefixes, or suffixes.

Once you have determined the fields that belong together, you can apply the Extract Class refactoring. If the component makes sense as a sub-class, Extract Subclass is also a candidate, and is often faster.

Many Parameters

Methods with many parameters are not only hard to understand, but their parameters also often become inconsistent when you need more, or different data.

There are several approaches to avoid long parameter lists:

1
# coding: utf-8
2
# Copyright (c) 2023 Stefan Bender
3
#
4
# This file is part of pyeppaurora.
5
# pyeppaurora is free software: you can redistribute it or modify
6
# it under the terms of the GNU General Public License as published
7
# by the Free Software Foundation, version 2.
8
# See accompanying LICENSE file or http://www.gnu.org/licenses/gpl-2.0.html.
9
"""Empirical model for auroral ionization rates
10
11
Implements the empirical model for auroral ionization,
12
derived from SSUSI UV observations [SB23]_.
13
14
.. [SB23] Bender et al., in prep., 2023, preprint: https://doi.org/10.48550/arXiv.2312.11130
15
"""
16
from logging import warning as warn
17
from os import path
18
from pkg_resources import resource_filename
19
20
import numpy as np
21
import xarray as xr
22
23
__all__ = [
24
	"ssusiq2023",
25
	"ssusiq2023_coeffs",
26
]
27
28
COEFF_FILE = "SSUSI_IRgrid_coeffs_f17f18.nc"
29
COEFF_PATH = resource_filename(__name__, path.join("data", COEFF_FILE))
30
31
32
def _interp(ds, method="linear", method_non_numeric="nearest", **kwargs):
33
	"""Fix `xarray` interpolation with non-numeric variables
34
	"""
35
	v_n = sorted(
36
		filter(lambda _v: np.issubdtype(ds[_v].dtype, np.number), ds)
37
	)
38
	v_nn = sorted(set(ds) - set(v_n))
39
	ds_n = ds[v_n].interp(method=method, **kwargs)
40
	ds_nn = ds[v_nn].sel(method=method_non_numeric, **kwargs)
41
	# override coordinates for `merge()`
42
	ds_nn = ds_nn.assign_coords(**ds_n.coords)
43
	return xr.merge([ds_n, ds_nn], join="left")
44
45
46
def ssusiq2023_coeffs(file=None):
47
	"""SSUSI ionization rate model coefficients
48
49
	Returns the fitted ionization rate model coefficents as
50
	read from the coefficient netcdf file.
51
52
	Parameters
53
	----------
54
	file: str, optional
55
		The coefficient file, defaults to the packaged coefficients.
56
57
	Returns
58
	-------
59
	coeffs: `xarray.Dataset`
60
		The default fitted model coefficients as read from the file.
61
	"""
62
	return xr.open_dataset(
63
		file or COEFF_PATH, decode_times=False, engine="h5netcdf"
64
	)
65
66
67
def ssusiq2023(
68
	gmlat,
69
	mlt,
70
	alt,
71
	sw_coeffs,
72
	coeff_ds=None,
73
	interpolate=False,
74
	method="linear",
75
	return_var=False,
76
):
77
	u"""
78
	Parameters
79
	----------
80
	gmlat: float
81
		Geomagnetic latitude in [degrees].
82
	mlt: float
83
		Magnetic local time in [hours].
84
	alt: float
85
		Altitude in [km]
86
	sw_coeffs: array_like or `xarray.DataArray`
87
		The space weather index values to use (for the requested time(s)),
88
		should be of shape (N, M) with N = number of proxies, currently 4:
89
		[Kp, PC, Ap, log(f10.7_81ctr_obs)].
90
		The `xarray.DataArray` should have a dimension named "proxy" with
91
		matching coordinates:
92
		["Kp", "PC", "Ap", "log_f107_81ctr_obs"]
93
		All the other dimensions will be broadcasted.
94
	coeff_ds: `xarray.Dataset`, optional (default: None)
95
		Dataset with the model coefficients, `None` uses the packaged version.
96
	interpolate: bool, optional (default: False)
97
		If `True`, uses bilinear interpolate in MLT and geomagnetic latitude,
98
		using periodic (24h) boundary conditions in MLT. Otherwise, the closest
99
		MLT/geomagnetic latitude bin will be selected.
100
	method: str, optional (default: "linear")
101
		Interpolation method to use, see `scipy.interpolate.interpn` for options.
102
		Only used if `interpolate` is `True`.
103
	return_var: bool, optional (default: False)
104
		If `True`, returns the predicted variance in addition to the values,
105
		otherwise only the mean prediction is returned.
106
107
	Returns
108
	-------
109
	q: `xarray.DataArray`
110
		log(q), where q is the ionization rate in [cm⁻³ s⁻¹]
111
		if `return_var` is False.
112
	q, var(q): tuple of `xarray.DataArray`s
113
		log(q) and var(log(q)) where q is the ionization rate in [cm⁻³ s⁻¹]
114
		if `return_var` is True.
115
	"""
116
	coeff_ds = coeff_ds or ssusiq2023_coeffs()
117
	coeff_sel = coeff_ds.sel(altitude=alt)
118
	if interpolate:
119
		_ds_m = coeff_sel.assign_coords(mlt=coeff_sel.mlt - 24)
120
		_ds_p = coeff_sel.assign_coords(mlt=coeff_sel.mlt + 24)
121
		_ds_mp = xr.concat([_ds_m, coeff_sel, _ds_p], dim="mlt")
122
		# square the standard deviation for interpolation
123
		_ds_mp["beta_var"] = _ds_mp["beta_std"]**2
124
		coeff_sel = _interp(
125
			_ds_mp,
126
			latitude=gmlat, mlt=mlt,
127
			method=method,
128
		)
129
		# and square root back to get the standard deviation
130
		coeff_sel["beta_std"] = np.sqrt(coeff_sel["beta_var"])
131
	else:
132
		coeff_sel = coeff_sel.sel(latitude=gmlat, mlt=mlt, method="nearest")
133
134
	# Determine if `xarray` read bytes or strings to
135
	# match the correct name in the proxy names.
136
	# Default is plain strings.
137
	offset = "offset"
138
	if isinstance(coeff_sel.proxy.values[0], bytes):
139
		offset = offset.encode()
140
	have_offset = offset in coeff_sel.proxy.values
141
142
	# prepare the coefficients (array) as a `xarray.DataArray`
143
	if isinstance(sw_coeffs, xr.DataArray):
144
		if have_offset:
145
			ones = xr.ones_like(sw_coeffs.isel(proxy=0))
146
			ones = ones.assign_coords(proxy="offset")
147
			sw_coeffs = xr.concat([sw_coeffs, ones], dim="proxy")
148
		sw_coeffs = sw_coeffs.sel(proxy=coeff_sel.proxy.astype(sw_coeffs.proxy.dtype))
149
	else:
150
		sw_coeffs = np.atleast_2d(sw_coeffs)
151
		if have_offset:
152
			aix = sw_coeffs.shape.index(len(coeff_sel.proxy.values) - 1)
153
			if aix != 0:
154
				warn(
155
					"Automatically changing axis. "
156
					"This is ambiguous, to remove the ambiguity, "
157
					"make sure that the different indexes (proxies) "
158
					"are ordered along the zero-th axis in multi-"
159
					"dimensional settings. I.e. each row corresponds "
160
					"to a different index, Kp, PC, Ap, etc."
161
				)
162
				sw_coeffs = sw_coeffs.swapaxes(aix, 0)
163
			sw_coeffs = np.vstack([sw_coeffs, np.ones(sw_coeffs.shape[1])])
164
		else:
165
			aix = sw_coeffs.shape.index(len(coeff_sel.proxy.values))
166
			if aix != 0:
167
				warn(
168
					"Automatically changing axis. "
169
					"This is ambiguous, to remove the ambiguity, "
170
					"make sure that the different indexes (proxies) "
171
					"are ordered along the zero-th axis in multi-"
172
					"dimensional settings. I.e. each row corresponds "
173
					"to a different index, Kp, PC, Ap, etc."
174
				)
175
				sw_coeffs = sw_coeffs.swapaxes(aix, 0)
176
		extra_dims = ["dim_{0}".format(_d) for _d in range(sw_coeffs.ndim - 1)]
177
		sw_coeffs = xr.DataArray(
178
			sw_coeffs,
179
			dims=["proxy"] + extra_dims,
180
			coords={"proxy": coeff_sel.proxy.values},
181
		)
182
183
	# Calculate model (mean) values from `beta`
184
	# fill NaNs with zero for `.dot()`
185
	coeffs = coeff_sel.beta.fillna(0.)
186
	q = coeffs.dot(sw_coeffs)
187
	q = q.rename("log_q")
188
	q.attrs = {
189
		"long_name": "natural logarithm of ionization rate",
190
		"units": "log(cm-3 s-1)",
191
	}
192
	if not return_var:
193
		return q
194
195
	# Calculate variance of the model from `beta_std`
196
	# fill NaNs with zero for `.dot()`
197
	coeffv = coeff_sel.beta_std.fillna(0.)**2
198
	q_var = coeffv.dot(sw_coeffs**2)
199
	if "sigma2" in coeff_sel.data_vars:
200
		# if available, add the posterior variance
201
		# to get the full posterior predictive variance
202
		q_var = coeff_sel["sigma2"] + q_var
203
	q_var = q_var.rename("var_log_q")
204
	q_var.attrs = {
205
		"long_name": "variance of the natural logarithm of ionization rate",
206
		"units": "1",
207
	}
208
	return q, q_var
209