Skip to content

Commit d6536c9

Browse files
adriessecwhanse
authored andcommitted
Implement IEC 61853 IAM calculations (#752)
* Insert three function placeholders in pvsystem.py * transform docstring from matlab * Code for iam_martin and iam_interp using greek symbol for theta * Minor improvements. * Tests for iam_martin_ruiz() and iam_interp() * Removed undesirables; polished docstrings and code here and there. * Final stickler. * Update v0.7.0.rst * Fix many oversights (also an old one). * Final(?) documentation touch-ups. * Fix tests and nan handling. * Final(?) test fixes. * Update v0.7.0.rst
1 parent 75d092e commit d6536c9

File tree

3 files changed

+246
-2
lines changed

3 files changed

+246
-2
lines changed

docs/sphinx/source/whatsnew/v0.7.0.rst

Lines changed: 7 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -10,10 +10,16 @@ compatibility notes.
1010
**Python 2.7 support ended on June 1, 2019.** (:issue:`501`)
1111
**Minimum numpy version is now 1.10.4. Minimum pandas version is now 0.18.1.**
1212

13+
Enhancements
14+
~~~~~~~~~~~~
15+
* Created two new incidence angle modifier functions: :py:func:`pvlib.pvsystem.iam_martin_ruiz`
16+
and :py:func:`pvlib.pvsystem.iam_interp`. (:issue:`751`)
17+
1318
Bug fixes
1419
~~~~~~~~~
1520
* Fix handling of keyword arguments in `forecasts.get_processed_data`.
1621
(:issue:`745`)
22+
* Fix output as Series feature in :py:func:`pvlib.pvsystem.ashraeiam`.
1723

1824
Testing
1925
~~~~~~~
@@ -26,4 +32,5 @@ Contributors
2632
* Mark Campanellli (:ghuser:`markcampanelli`)
2733
* Will Holmgren (:ghuser:`wholmgren`)
2834
* Oscar Dowson (:ghuser:`odow`)
35+
* Anton Driesse (:ghuser:`adriesse`)
2936
* Alexander Morgan (:ghuser:`alexandermorgan`)

pvlib/pvsystem.py

Lines changed: 157 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -17,7 +17,7 @@
1717
import pandas as pd
1818

1919
from pvlib import atmosphere, irradiance, tools, singlediode as _singlediode
20-
from pvlib.tools import _build_kwargs
20+
from pvlib.tools import _build_kwargs, cosd
2121
from pvlib.location import Location
2222

2323

@@ -946,7 +946,7 @@ def ashraeiam(aoi, b=0.05):
946946
iam = np.where(aoi_gte_90, 0, iam)
947947
iam = np.maximum(0, iam)
948948

949-
if isinstance(iam, pd.Series):
949+
if isinstance(aoi, pd.Series):
950950
iam = pd.Series(iam, index=aoi.index)
951951

952952
return iam
@@ -1064,6 +1064,161 @@ def physicaliam(aoi, n=1.526, K=4., L=0.002):
10641064
return iam
10651065

10661066

1067+
def iam_martin_ruiz(aoi, a_r=0.16):
1068+
'''
1069+
Determine the incidence angle modifier (iam) using the Martin
1070+
and Ruiz incident angle model.
1071+
1072+
Parameters
1073+
----------
1074+
aoi : numeric, degrees
1075+
The angle of incidence between the module normal vector and the
1076+
sun-beam vector in degrees. Theta must be a numeric scalar or vector.
1077+
iam is 0 where |aoi| > 90.
1078+
1079+
a_r : numeric
1080+
The angular losses coefficient described in equation 3 of [1].
1081+
This is an empirical dimensionless parameter. Values of a_r are
1082+
generally on the order of 0.08 to 0.25 for flat-plate PV modules.
1083+
a_r must be a positive numeric scalar or vector (same length as aoi).
1084+
1085+
Returns
1086+
-------
1087+
iam : numeric
1088+
The incident angle modifier(s)
1089+
1090+
Notes
1091+
-----
1092+
iam_martin_ruiz calculates the incidence angle modifier (iam)
1093+
as described by Martin and Ruiz in [1]. The information
1094+
required is the incident angle (aoi) and the angular losses
1095+
coefficient (a_r). Please note that [1] has a corrigendum which makes
1096+
the document much simpler to understand.
1097+
1098+
The incident angle modifier is defined as
1099+
[1-exp(-cos(aoi/ar))] / [1-exp(-1/ar)], which is
1100+
presented as AL(alpha) = 1 - IAM in equation 4 of [1]. Thus IAM is
1101+
equal to 1 at aoi = 0, and equal to 0 at aoi = 90. This equation is only
1102+
valid for -90 <= aoi <= 90, therefore iam must be constrained to 0.0
1103+
beyond this range.
1104+
1105+
References
1106+
----------
1107+
[1] N. Martin and J. M. Ruiz, "Calculation of the PV modules angular
1108+
losses under field conditions by means of an analytical model", Solar
1109+
Energy Materials & Solar Cells, vol. 70, pp. 25-38, 2001.
1110+
1111+
[2] N. Martin and J. M. Ruiz, "Corrigendum to 'Calculation of the PV
1112+
modules angular losses under field conditions by means of an
1113+
analytical model'", Solar Energy Materials & Solar Cells, vol. 110,
1114+
pp. 154, 2013.
1115+
1116+
See Also
1117+
--------
1118+
physicaliam
1119+
ashraeiam
1120+
iam_interp
1121+
'''
1122+
# Contributed by Anton Driesse (@adriesse), PV Performance Labs. July, 2019
1123+
1124+
aoi_input = aoi
1125+
1126+
aoi = np.asanyarray(aoi)
1127+
a_r = np.asanyarray(a_r)
1128+
1129+
if np.any(np.less_equal(a_r, 0)):
1130+
raise RuntimeError("The parameter 'a_r' cannot be zero or negative.")
1131+
1132+
with np.errstate(invalid='ignore'):
1133+
iam = (1 - np.exp(-cosd(aoi) / a_r)) / (1 - np.exp(-1 / a_r))
1134+
iam = np.where(np.abs(aoi) >= 90.0, 0.0, iam)
1135+
1136+
if isinstance(aoi_input, pd.Series):
1137+
iam = pd.Series(iam, index=aoi_input.index)
1138+
1139+
return iam
1140+
1141+
1142+
def iam_interp(aoi, theta_ref, iam_ref, method='linear', normalize=True):
1143+
'''
1144+
Determine the incidence angle modifier (iam) by interpolating a set of
1145+
reference values, which are usually measured values.
1146+
1147+
Parameters
1148+
----------
1149+
aoi : numeric, degrees
1150+
The angle of incidence between the module normal vector and the
1151+
sun-beam vector in degrees.
1152+
1153+
theta_ref : numeric, degrees
1154+
Vector of angles at which the iam is known.
1155+
1156+
iam_ref : numeric, unitless
1157+
iam values for each angle in theta_ref.
1158+
1159+
method : str, default 'linear'
1160+
Specifies the interpolation method.
1161+
Useful options are: 'linear', 'quadratic','cubic'.
1162+
See scipy.interpolate.interp1d for more options.
1163+
1164+
normalize : boolean
1165+
When true, the interpolated values are divided by the interpolated
1166+
value at zero degrees. This ensures that the iam at normal
1167+
incidence is equal to 1.0.
1168+
1169+
Returns
1170+
-------
1171+
iam : numeric
1172+
The incident angle modifier(s)
1173+
1174+
Notes:
1175+
------
1176+
theta_ref must have two or more points and may span any range of angles.
1177+
Typically there will be a dozen or more points in the range 0-90 degrees.
1178+
iam beyond the range of theta_ref are extrapolated, but constrained to be
1179+
non-negative.
1180+
1181+
The sign of aoi is ignored; only the magnitude is used.
1182+
1183+
See Also
1184+
--------
1185+
physicaliam
1186+
ashraeiam
1187+
iam_martin_ruiz
1188+
'''
1189+
# Contributed by Anton Driesse (@adriesse), PV Performance Labs. July, 2019
1190+
1191+
from scipy.interpolate import interp1d
1192+
1193+
# Scipy doesn't give the clearest feedback, so check number of points here.
1194+
MIN_REF_VALS = {'linear': 2, 'quadratic': 3, 'cubic': 4, 1: 2, 2: 3, 3: 4}
1195+
1196+
if len(theta_ref) < MIN_REF_VALS.get(method, 2):
1197+
raise ValueError("Too few reference points defined "
1198+
"for interpolation method '%s'." % method)
1199+
1200+
if np.any(np.less(iam_ref, 0)):
1201+
raise ValueError("Negative value(s) found in 'iam_ref'. "
1202+
"This is not physically possible.")
1203+
1204+
interpolator = interp1d(theta_ref, iam_ref, kind=method,
1205+
fill_value='extrapolate')
1206+
aoi_input = aoi
1207+
1208+
aoi = np.asanyarray(aoi)
1209+
aoi = np.abs(aoi)
1210+
iam = interpolator(aoi)
1211+
iam = np.clip(iam, 0, None)
1212+
1213+
if normalize:
1214+
iam /= interpolator(0)
1215+
1216+
if isinstance(aoi_input, pd.Series):
1217+
iam = pd.Series(iam, index=aoi_input.index)
1218+
1219+
return iam
1220+
1221+
10671222
def calcparams_desoto(effective_irradiance, temp_cell,
10681223
alpha_sc, a_ref, I_L_ref, I_o_ref, R_sh_ref, R_s,
10691224
EgRef=1.121, dEgdT=-0.0002677,

pvlib/test/test_pvsystem.py

Lines changed: 82 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -148,6 +148,88 @@ def test_PVSystem_physicaliam(mocker):
148148
assert iam < 1.
149149

150150

151+
def test_iam_martin_ruiz():
152+
153+
aoi = 45.
154+
a_r = 0.16
155+
expected = 0.98986965
156+
157+
# will fail of default values change
158+
iam = pvsystem.iam_martin_ruiz(aoi)
159+
assert_allclose(iam, expected)
160+
# will fail of parameter names change
161+
iam = pvsystem.iam_martin_ruiz(aoi=aoi, a_r=a_r)
162+
assert_allclose(iam, expected)
163+
164+
a_r = 0.18
165+
aoi = [-100, -60, 0, 60, 100, np.nan, np.inf]
166+
expected = [0.0, 0.9414631, 1.0, 0.9414631, 0.0, np.nan, 0.0]
167+
168+
# check out of range of inputs as list
169+
iam = pvsystem.iam_martin_ruiz(aoi, a_r)
170+
assert_allclose(iam, expected, equal_nan=True)
171+
172+
# check out of range of inputs as array
173+
iam = pvsystem.iam_martin_ruiz(np.array(aoi), a_r)
174+
assert_allclose(iam, expected, equal_nan=True)
175+
176+
# check out of range of inputs as Series
177+
aoi = pd.Series(aoi)
178+
expected = pd.Series(expected)
179+
iam = pvsystem.iam_martin_ruiz(aoi, a_r)
180+
assert_series_equal(iam, expected)
181+
182+
# check exception clause
183+
with pytest.raises(RuntimeError):
184+
pvsystem.iam_martin_ruiz(0.0, a_r=0.0)
185+
186+
187+
@requires_scipy
188+
def test_iam_interp():
189+
190+
aoi_meas = [0.0, 45.0, 65.0, 75.0]
191+
iam_meas = [1.0, 0.9, 0.8, 0.6]
192+
193+
# simple default linear method
194+
aoi = 55.0
195+
expected = 0.85
196+
iam = pvsystem.iam_interp(aoi, aoi_meas, iam_meas)
197+
assert_allclose(iam, expected)
198+
199+
# simple non-default method
200+
aoi = 55.0
201+
expected = 0.8878062
202+
iam = pvsystem.iam_interp(aoi, aoi_meas, iam_meas, method='cubic')
203+
assert_allclose(iam, expected)
204+
205+
# check with all reference values
206+
aoi = aoi_meas
207+
expected = iam_meas
208+
iam = pvsystem.iam_interp(aoi, aoi_meas, iam_meas)
209+
assert_allclose(iam, expected)
210+
211+
# check normalization and Series
212+
aoi = pd.Series(aoi)
213+
expected = pd.Series(expected)
214+
iam_mult = np.multiply(0.9, iam_meas)
215+
iam = pvsystem.iam_interp(aoi, aoi_meas, iam_mult, normalize=True)
216+
assert_series_equal(iam, expected)
217+
218+
# check beyond reference values
219+
aoi = [-45, 0, 45, 85, 90, 95, 100, 105, 110]
220+
expected = [0.9, 1.0, 0.9, 0.4, 0.3, 0.2, 0.1, 0.0, 0.0]
221+
iam = pvsystem.iam_interp(aoi, aoi_meas, iam_meas)
222+
assert_allclose(iam, expected)
223+
224+
# check exception clause
225+
with pytest.raises(ValueError):
226+
pvsystem.iam_interp(0.0, [0], [1])
227+
228+
# check exception clause
229+
with pytest.raises(ValueError):
230+
pvsystem.iam_interp(0.0, [0, 90], [1, -1])
231+
232+
151233
# if this completes successfully we'll be able to do more tests below.
152234
@pytest.fixture(scope="session")
153235
def sam_data():

0 commit comments

Comments
 (0)