-
Notifications
You must be signed in to change notification settings - Fork 1.1k
coefficient estimation method following DeSoto(2006) #784
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
Changes from 34 commits
c5287ea
5a39fea
7b24201
c16f1d0
e7ef571
cdb2a73
0cfe13f
0c27829
d8653cd
bbb4580
6b8fadf
121cb06
47e1be4
c1b62c4
4e62d1a
747ed48
9ea2ed8
bd39dbd
ae5c8be
2d4f8f3
485dae6
77d88a2
6fcafc7
2d53d5d
fdbf5ec
0e1b4b3
d2d8c45
a70debb
e01f262
17f8617
e14ad40
32049b0
bfe3994
7ff8a96
7f08d80
27ab961
79ef920
5e9cf8d
c02e74f
fa41eb8
753d312
b1b405c
22d613c
ae19a68
6913a52
05051e1
d4772c6
d821494
9d23552
5037da6
File filter
Filter by extension
Conversations
Jump to
Diff view
Diff view
There are no files selected for viewing
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -1,4 +1,4 @@ | ||
.. _whatsnew_0700: | ||
.. _whatsnew_0700: | ||
|
||
v0.7.0 (MONTH DAY, YEAR) | ||
------------------------ | ||
|
@@ -96,6 +96,8 @@ Enhancements | |
the single diode equation to an IV curve. | ||
* Add :py:func:`~pvlib.ivtools.fit_sdm_cec_sam`, a wrapper for the CEC single | ||
diode model fitting function '6parsolve' from NREL's System Advisor Model. | ||
* Add :py:func:`~pvlib.ivtools.fit_sdm_desoto`, a method to fit the single | ||
cwhanse marked this conversation as resolved.
Show resolved
Hide resolved
cwhanse marked this conversation as resolved.
Show resolved
Hide resolved
|
||
diode equation to the typical specifications given in manufacturers datasheets. | ||
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. Probably should say “DeSoto single diode model” instead of “single diode equation”. |
||
|
||
Bug fixes | ||
~~~~~~~~~ | ||
|
Original file line number | Diff line number | Diff line change |
---|---|---|
|
@@ -6,6 +6,8 @@ | |
""" | ||
|
||
import numpy as np | ||
from scipy.optimize import root | ||
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. The |
||
from scipy import constants | ||
|
||
|
||
def fit_sdm_cec_sam(celltype, v_mp, i_mp, v_oc, i_sc, alpha_sc, beta_voc, | ||
|
@@ -262,6 +264,163 @@ def fit_sde_sandia(voltage, current, v_oc=None, i_sc=None, v_mp_i_mp=None, | |
v_oc) | ||
|
||
|
||
def fit_sdm_desoto(v_mp, i_mp, v_oc, i_sc, alpha_sc, beta_voc, | ||
cells_in_series, EgRef=1.121, dEgdT=-0.0002677, | ||
temp_ref=25, irrad_ref=1000, solver_method='lm'): | ||
""" | ||
Calculates the parameters for the De Soto single diode model using the | ||
procedure described in [1]. This procedure has the | ||
advantage of using common specifications given by manufacturers in the | ||
datasheets of PV modules. | ||
The parameters returned by this function can be used by | ||
pvsystem.calcparams_desoto to calculate the values at different | ||
irradiance and cell temperature. | ||
|
||
Parameters | ||
---------- | ||
v_mp: float | ||
cwhanse marked this conversation as resolved.
Show resolved
Hide resolved
|
||
Module voltage at the maximum-power point at reference conditions [V]. | ||
i_mp: float | ||
Module current at the maximum-power point at reference conditions [A]. | ||
v_oc: float | ||
Open-circuit voltage at reference conditions [V]. | ||
i_sc: float | ||
Short-circuit current at reference conditions [A]. | ||
alpha_sc: float | ||
The short-circuit current (i_sc) temperature coefficient of the | ||
module [A/K]. | ||
beta_voc: float | ||
The open-circuit voltage (v_oc) temperature coefficient of the | ||
module [V/K]. | ||
cells_in_series: float | ||
cwhanse marked this conversation as resolved.
Show resolved
Hide resolved
|
||
Number of cell in the module. | ||
tylunel marked this conversation as resolved.
Show resolved
Hide resolved
|
||
EgRef: float, default 1.121 eV - value for silicon | ||
Energy of bandgap of semi-conductor used [eV] | ||
dEgdT: float, default -0.0002677 - value for silicon | ||
Variation of bandgap according to temperature [eV/K] | ||
temp_ref: float, default 25 | ||
Reference temperature condition [C] | ||
irrad_ref: float, default 1000 | ||
Reference irradiance condition [W/m2] | ||
solver_method: str, default 'lm' | ||
See help of scipy.optimize.root() for more details. | ||
tylunel marked this conversation as resolved.
Show resolved
Hide resolved
|
||
|
||
Returns | ||
------- | ||
Dictionary with the following elements: | ||
|
||
* I_L_ref: float | ||
Light-generated current at reference conditions [A] | ||
* I_o_ref: float | ||
Diode saturation current at reference conditions [A] | ||
* R_s: float | ||
Series resistance [ohms] | ||
Note that '_ref' is not mentionned in the name because | ||
cwhanse marked this conversation as resolved.
Show resolved
Hide resolved
|
||
this resistance is not sensible to the conditions of test. | ||
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. I don't readily understand what this note means. I always thought the missing There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more.
|
||
* R_sh_ref: float | ||
Shunt resistance at reference conditions [ohms]. | ||
* a_ref: float | ||
Modified ideality factor at reference conditions. | ||
The product of the usual diode ideality factor (n, unitless), | ||
number of cells in series (Ns), and cell thermal voltage at | ||
specified effective irradiance and cell temperature. | ||
* alpha_sc: float | ||
Caution!: Different from the input because of the unit. | ||
The short-circuit current (i_sc) temperature coefficient of the | ||
module [A/K]. | ||
* EgRef: float | ||
Energy of bandgap of semi-conductor used [eV] | ||
* dEgdT: float | ||
Variation of bandgap according to temperature [eV/K] | ||
* irrad_ref: float | ||
Reference irradiance condition [W/m2] | ||
* temp_ref: float | ||
Reference temperature condition [C] | ||
|
||
References | ||
---------- | ||
[1] W. De Soto et al., "Improvement and validation of a model for | ||
photovoltaic array performance", Solar Energy, vol 80, pp. 78-88, | ||
2006. | ||
|
||
[2] John A Duffie ,William A Beckman, "Solar Engineering of Thermal | ||
cwhanse marked this conversation as resolved.
Show resolved
Hide resolved
|
||
Processes", Wiley, 2013 | ||
""" | ||
# Constants | ||
k = constants.value('Boltzmann constant in eV/K') | ||
Tref = temp_ref + 273.15 # [K] | ||
|
||
def pv_fct(params, specs): | ||
cwhanse marked this conversation as resolved.
Show resolved
Hide resolved
|
||
"""Evaluates the systems of equations used to solve for the single | ||
diode equation parameters. | ||
To avoid the confusion in names with variables of container | ||
cwhanse marked this conversation as resolved.
Show resolved
Hide resolved
|
||
function the '_' of the variables were removed. | ||
""" | ||
# six input known variables | ||
Isc, Voc, Imp, Vmp, betaoc, alphasc = specs | ||
|
||
# five parameters vector to find | ||
IL, Io, a, Rsh, Rs = params | ||
|
||
# five equation vector | ||
y = [0, 0, 0, 0, 0] | ||
|
||
# 1st equation - short-circuit - eq(3) in [1] | ||
y[0] = Isc - IL + Io*np.expm1(Isc*Rs/a) + Isc*Rs/Rsh | ||
# 2nd equation - open-circuit Tref - eq(4) in [1] | ||
y[1] = -IL + Io*np.expm1(Voc/a) + Voc/Rsh | ||
# 3rd equation - Imp & Vmp - eq(5) in [1] | ||
y[2] = Imp - IL + Io*np.expm1((Vmp+Imp*Rs)/a) + \ | ||
(Vmp+Imp*Rs)/Rsh | ||
# 4th equation - Pmp derivated=0 - eq23.2.6 in [2] | ||
# caution: eq(6) in [1] has a sign error | ||
y[3] = Imp - Vmp * ((Io/a)*np.exp((Vmp+Imp*Rs)/a) + 1.0/Rsh) / \ | ||
(1.0 + (Io*Rs/a)*np.exp((Vmp+Imp*Rs)/a) + Rs/Rsh) | ||
# 5th equation - open-circuit T2 - eq (4) at temperature T2 in [1] | ||
T2 = Tref + 2 | ||
Voc2 = (T2 - Tref)*betaoc + Voc # eq (7) in [1] | ||
a2 = a*T2/Tref # eq (8) in [1] | ||
IL2 = IL + alphasc*(T2-Tref) # eq (11) in [1] | ||
Eg2 = EgRef*(1 + dEgdT*(T2-Tref)) # eq (10) in [1] | ||
Io2 = Io * (T2/Tref)**3 * np.exp(1/k * (EgRef/Tref-Eg2/T2)) # eq (9) | ||
y[4] = -IL2 + Io2*np.expm1(Voc2/a2) + Voc2/Rsh # eq (4) at T2 | ||
|
||
return y | ||
|
||
# initial guesses of variables for computing convergence: | ||
# Values are taken from [2], p753 | ||
Rsh_i = 100.0 | ||
a_i = 1.5*k*Tref*cells_in_series | ||
IL_i = i_sc | ||
Io_i = i_sc * np.exp(-v_oc/a_i) | ||
Rs_i = (a_i*np.log1p((IL_i-i_mp)/Io_i) - v_mp)/i_mp | ||
# params_i : initial values vector | ||
params_i = np.array([IL_i, Io_i, a_i, Rsh_i, Rs_i]) | ||
cwhanse marked this conversation as resolved.
Show resolved
Hide resolved
|
||
|
||
# specs of module | ||
specs = np.array([i_sc, v_oc, i_mp, v_mp, beta_voc, alpha_sc]) | ||
cwhanse marked this conversation as resolved.
Show resolved
Hide resolved
|
||
|
||
# computing | ||
result = root(pv_fct, x0=params_i, args=specs, method=solver_method) | ||
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. Have you attempted to compute the Jacobian of This might be both faster and more robust, but it also can easily be left as a followup improvement. There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. Not necessarily true. I believe the secant method is more stable (*) than the Newton-Raphson, and even though it converges slower, it might still be faster, if the Jacobian is computationally expensive to calculate. (*) I believe secant or Broyden's is less sensitive to poor initial guess, discontinuity, and relatively flat or stiff systems, which can cause N-R to get stuck in an infinite loop, but I can't find a good reference for this, so maybe I'm wrong, sorry There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. BTW, if you do need the derivatives, check There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more.
Agree, let's tackle that with a follow-on PR. |
||
|
||
if result.success: | ||
sdm_params = result.x | ||
else: | ||
raise RuntimeError('Parameter estimation failed') | ||
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. Maybe include the full There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. Since we're not exposing most of the scipy.optimize.root input options, I'm not in favor of returning the |
||
|
||
# results | ||
return {'I_L_ref': sdm_params[0], | ||
'I_o_ref': sdm_params[1], | ||
'a_ref': sdm_params[2], | ||
'R_sh_ref': sdm_params[3], | ||
'R_s': sdm_params[4], | ||
'alpha_sc': alpha_sc, | ||
'EgRef': EgRef, | ||
'dEgdT': dEgdT, | ||
'irrad_ref': irrad_ref, | ||
'temp_ref': temp_ref} | ||
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. I personally like to be able to inspect the full There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. Ok. There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. I don't think it is worth testing the content of However, your question led me to look at the default settings for |
||
|
||
|
||
def _find_mp(voltage, current): | ||
""" | ||
Finds voltage and current at maximum power point. | ||
|
Original file line number | Diff line number | Diff line change |
---|---|---|
|
@@ -28,6 +28,13 @@ def get_cec_params_cansol_cs5p_220p(): | |
'R_sh_ref': 837.51, 'R_s': 1.004, 'Adjust': 2.3}} | ||
|
||
|
||
@pytest.fixture | ||
def get_test_specs_params(): | ||
cwhanse marked this conversation as resolved.
Show resolved
Hide resolved
|
||
"""Specifications of module Kyocera KU270-6MCA""" | ||
return {'v_mp': 31.0, 'i_mp': 8.71, 'v_oc': 38.3, | ||
'i_sc': 9.43, 'alpha_sc': 0.005658, 'beta_voc': -0.13788} | ||
|
||
|
||
@requires_scipy | ||
def test_fit_sde_sandia(get_test_iv_params, get_bad_iv_curves): | ||
test_params = get_test_iv_params | ||
|
@@ -102,6 +109,29 @@ def test_fit_sdm_cec_sam(get_cec_params_cansol_cs5p_220p): | |
cells_in_series=1, temp_ref=25) | ||
|
||
|
||
@requires_scipy | ||
def test_fit_sdm_desoto(get_test_specs_params): | ||
result = ivtools.fit_sdm_desoto(cells_in_series=60, | ||
tylunel marked this conversation as resolved.
Show resolved
Hide resolved
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. move There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. cells_in_series is not in the get_test_specs_params fixture because it changes at line 129 of the test in order to cause the RuntimeError. I feel like it would be more complicated to add it the fixture and to change it later to have the RuntimeError. What do you think ? There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. I'd prefer to see it added to the fixture and then modified in the function that causes the |
||
**get_test_specs_params) | ||
result_expected = {'I_L_ref': 9.452324509050774, | ||
'I_o_ref': 3.2246097466679494e-10, | ||
'a_ref': 1.5912875522463978, | ||
'R_sh_ref': 125.79869968976132, | ||
'R_s': 0.2978148476600744, | ||
'alpha_sc': 0.005658, | ||
'EgRef': 1.121, | ||
'dEgdT': -0.0002677, | ||
'irrad_ref': 1000, | ||
'temp_ref': 25} | ||
assert np.allclose(pd.Series(result), pd.Series(result_expected), | ||
rtol=1e-4) | ||
cwhanse marked this conversation as resolved.
Show resolved
Hide resolved
|
||
with pytest.raises(RuntimeError): | ||
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. Could check that the exception raised is the exception we're after. Could do it like this (learned from @wholmgren quite recently):
|
||
ivtools.fit_sdm_desoto(cells_in_series=10, | ||
v_mp=31.0, i_mp=8.71, v_oc=38.3, | ||
i_sc=9.43, alpha_sc=0.06, | ||
beta_voc=-0.36) | ||
|
||
|
||
@pytest.fixture | ||
def get_bad_iv_curves(): | ||
# v1, i1 produces a bad value for I0_voc | ||
|
Uh oh!
There was an error while loading. Please reload this page.