Skip to content

Commit 10d54de

Browse files
tylunelcwhanse
authored andcommitted
coefficient estimation method following DeSoto(2006) (#784)
* function getparams_desoto added to pvsystem module. This commit is just the copy from my previous work, slighty modified for removing the PEP8 warnings. * getparams_desoto moved from pvsystem to singlediode. - getparams_desoto renamed in getparams_from_specs - some initializing values were changed in accordance with Duffie&Beckman2013 * - Modification of getparams_from_specs so as it follows better the procedure of DeSoto&al(2006). * - Bug corrected in DeSoto(2006) procedure - test_singlediode completed with a beginning of test_getparams_from_specs * - test_getparams_from_specs_desoto() finished - function renamed as above * Not sure of all changes brought by this commit because of holidays. - function 'from_epw' added in location - 1 modification in singlediode - 1 modification in test_singlediode * - function '_parse_raw_sam_df' modified. The parser engine is now defined on 'python'. If not the pd.read_csv cannot work with me. * - ModelChain attribute 'diode_params' transformed from tuple containing pd.Series to DataFrame. Makes the use of diode_params easier for further calculations. * - singlediode.get_params_from_specs_desoto() output changed. 'a_ref' is changed by 'nNsVth_ref' * read_epw changed. A line has been added to convert the precipitable water from mm to cm, in order to be compatible with other functions of pvlib * - read_epw changed. If condition added to make the conversion only in the case of TMY3 * - argument diode_params changed from tuple to pd.DataFrame * - get_params_from_specs_desoto removed from singlediode.py * - function fit_sdm_desoto added. Still need to be formatted * - change on type of self.diode_params removed. Go check on branch diode_params_in_df for seeing it * - Function 'fit_sdm_desoto' cleaned and variables names named as in 'fit_sdm_cec_sam' * - all changes made on other files than ivtools.py removed (cleaning for comparing before PR) * - other differences cleaned * - renaming of one variable and minor documentation modifications * - Beginning of writting of test_fit_sdm_desoto. Coverage around 90-95% I think * -minor format changes * - changes made according to feedbacks of markcampanelli * - some cleaning on fit_sdm_desoto to make it more readable - tests completed * - minor code cleaning * - check on importation of scipy removed * - minor cleaning * - attempt to reach 100% coverage * - description added in docs/sphinx/source/whatsnew/v0.7.0.rst * - changes made according to cwhanse review. Except alpha_sc and beta_voc still in %/K * - minor correction and adaptation of test * - sign correction on 3rd equation * - changes on units of alpha_sc and beta_voc inputs. Now in A/K and V/K rather than %/K. - 'celltype' input replaced by EgRef and dEgdT, with values of Si as default * - other line of test added for better coverage * - some changes to try feedbacks of cwhanse and markcampanelli, not finished * -OptimizeResult added in output - solver_method removed - docstring modified - result.message included in message raised by RuntimeError * - includes all feedbacks made on the 21/10, except moving of pv_fct() to the module level * - pv_fct moved out of the fit_sdm_desoto function and renamed in _system_of_equations * - minor modification: Boltzmann k given in specs to avoid import of scipy in _system_of_equations * - cleaning and minor modifications to docstring * - references added to docstring in _system_of_equations * - modification for removing last lint error * - modification for removing last lint error * - modification for removing lint error * - integration of adriesse suggestions * - adding of tylunel to the list of contributors * - adding of mark requires_scipy
1 parent b99ded7 commit 10d54de

File tree

4 files changed

+241
-1
lines changed

4 files changed

+241
-1
lines changed

docs/sphinx/source/api.rst

+1
Original file line numberDiff line numberDiff line change
@@ -303,6 +303,7 @@ Functions for fitting PV models
303303

304304
ivtools.fit_sde_sandia
305305
ivtools.fit_sdm_cec_sam
306+
ivtools.fit_sdm_desoto
306307

307308
Other
308309
-----

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

+4-1
Original file line numberDiff line numberDiff line change
@@ -1,4 +1,4 @@
1-
.. _whatsnew_0700:
1+
.. _whatsnew_0700:
22

33
v0.7.0 (MONTH DAY, YEAR)
44
------------------------
@@ -111,6 +111,8 @@ Enhancements
111111
the single diode equation to an IV curve.
112112
* Add :py:func:`~pvlib.ivtools.fit_sdm_cec_sam`, a wrapper for the CEC single
113113
diode model fitting function '6parsolve' from NREL's System Advisor Model.
114+
* Add :py:func:`~pvlib.ivtools.fit_sdm_desoto`, a method to fit the De Soto single
115+
diode model to the typical specifications given in manufacturers datasheets.
114116
* Add `timeout` to :py:func:`pvlib.iotools.get_psm3`.
115117

116118
Bug fixes
@@ -162,4 +164,5 @@ Contributors
162164
* Anton Driesse (:ghuser:`adriesse`)
163165
* Alexander Morgan (:ghuser:`alexandermorgan`)
164166
* Miguel Sánchez de León Peque (:ghuser:`Peque`)
167+
* Tanguy Lunel (:ghuser:`tylunel`)
165168
* Veronica Guo (:ghuser:`veronicaguo`)

pvlib/ivtools.py

+207
Original file line numberDiff line numberDiff line change
@@ -262,6 +262,147 @@ def fit_sde_sandia(voltage, current, v_oc=None, i_sc=None, v_mp_i_mp=None,
262262
v_oc)
263263

264264

265+
def fit_sdm_desoto(v_mp, i_mp, v_oc, i_sc, alpha_sc, beta_voc,
266+
cells_in_series, EgRef=1.121, dEgdT=-0.0002677,
267+
temp_ref=25, irrad_ref=1000, root_kwargs={}):
268+
"""
269+
Calculates the parameters for the De Soto single diode model using the
270+
procedure described in [1]. This procedure has the advantage of
271+
using common specifications given by manufacturers in the
272+
datasheets of PV modules.
273+
274+
The solution is found using the scipy.optimize.root() function,
275+
with the corresponding default solver method 'hybr'.
276+
No restriction is put on the fit variables, i.e. series
277+
or shunt resistance could go negative. Nevertheless, if it happens,
278+
check carefully the inputs and their units; alpha_sc and beta_voc are
279+
often given in %/K in manufacturers datasheets and should be given
280+
in A/K and V/K here.
281+
282+
The parameters returned by this function can be used by
283+
pvsystem.calcparams_desoto to calculate the values at different
284+
irradiance and cell temperature.
285+
286+
Parameters
287+
----------
288+
v_mp: float
289+
Module voltage at the maximum-power point at reference conditions [V].
290+
i_mp: float
291+
Module current at the maximum-power point at reference conditions [A].
292+
v_oc: float
293+
Open-circuit voltage at reference conditions [V].
294+
i_sc: float
295+
Short-circuit current at reference conditions [A].
296+
alpha_sc: float
297+
The short-circuit current (i_sc) temperature coefficient of the
298+
module [A/K].
299+
beta_voc: float
300+
The open-circuit voltage (v_oc) temperature coefficient of the
301+
module [V/K].
302+
cells_in_series: integer
303+
Number of cell in the module.
304+
EgRef: float, default 1.121 eV - value for silicon
305+
Energy of bandgap of semi-conductor used [eV]
306+
dEgdT: float, default -0.0002677 - value for silicon
307+
Variation of bandgap according to temperature [eV/K]
308+
temp_ref: float, default 25
309+
Reference temperature condition [C]
310+
irrad_ref: float, default 1000
311+
Reference irradiance condition [W/m2]
312+
root_kwargs: dictionary, default None
313+
Dictionary of arguments to pass onto scipy.optimize.root()
314+
315+
Returns
316+
-------
317+
Tuple of the following elements:
318+
319+
* Dictionary with the following elements:
320+
I_L_ref: float
321+
Light-generated current at reference conditions [A]
322+
I_o_ref: float
323+
Diode saturation current at reference conditions [A]
324+
R_s: float
325+
Series resistance [ohms]
326+
R_sh_ref: float
327+
Shunt resistance at reference conditions [ohms].
328+
a_ref: float
329+
Modified ideality factor at reference conditions.
330+
The product of the usual diode ideality factor (n, unitless),
331+
number of cells in series (Ns), and cell thermal voltage at
332+
specified effective irradiance and cell temperature.
333+
alpha_sc: float
334+
The short-circuit current (i_sc) temperature coefficient of the
335+
module [A/K].
336+
EgRef: float
337+
Energy of bandgap of semi-conductor used [eV]
338+
dEgdT: float
339+
Variation of bandgap according to temperature [eV/K]
340+
irrad_ref: float
341+
Reference irradiance condition [W/m2]
342+
temp_ref: float
343+
Reference temperature condition [C]
344+
* scipy.optimize.OptimizeResult
345+
Optimization result of scipy.optimize.root().
346+
See scipy.optimize.OptimizeResult for more details.
347+
348+
References
349+
----------
350+
[1] W. De Soto et al., "Improvement and validation of a model for
351+
photovoltaic array performance", Solar Energy, vol 80, pp. 78-88,
352+
2006.
353+
354+
[2] John A Duffie, William A Beckman, "Solar Engineering of Thermal
355+
Processes", Wiley, 2013
356+
"""
357+
358+
try:
359+
from scipy.optimize import root
360+
from scipy import constants
361+
except ImportError:
362+
raise ImportError("The fit_sdm_desoto function requires scipy.")
363+
364+
# Constants
365+
k = constants.value('Boltzmann constant in eV/K')
366+
Tref = temp_ref + 273.15 # [K]
367+
368+
# initial guesses of variables for computing convergence:
369+
# Values are taken from [2], p753
370+
Rsh_0 = 100.0
371+
a_0 = 1.5*k*Tref*cells_in_series
372+
IL_0 = i_sc
373+
Io_0 = i_sc * np.exp(-v_oc/a_0)
374+
Rs_0 = (a_0*np.log1p((IL_0-i_mp)/Io_0) - v_mp)/i_mp
375+
# params_i : initial values vector
376+
params_i = np.array([IL_0, Io_0, a_0, Rsh_0, Rs_0])
377+
378+
# specs of module
379+
specs = (i_sc, v_oc, i_mp, v_mp, beta_voc, alpha_sc, EgRef, dEgdT,
380+
Tref, k)
381+
382+
# computing with system of equations described in [1]
383+
optimize_result = root(_system_of_equations_desoto, x0=params_i,
384+
args=(specs,), **root_kwargs)
385+
386+
if optimize_result.success:
387+
sdm_params = optimize_result.x
388+
else:
389+
raise RuntimeError(
390+
'Parameter estimation failed:\n' + optimize_result.message)
391+
392+
# results
393+
return ({'I_L_ref': sdm_params[0],
394+
'I_o_ref': sdm_params[1],
395+
'a_ref': sdm_params[2],
396+
'R_sh_ref': sdm_params[3],
397+
'R_s': sdm_params[4],
398+
'alpha_sc': alpha_sc,
399+
'EgRef': EgRef,
400+
'dEgdT': dEgdT,
401+
'irrad_ref': irrad_ref,
402+
'temp_ref': temp_ref},
403+
optimize_result)
404+
405+
265406
def _find_mp(voltage, current):
266407
"""
267408
Finds voltage and current at maximum power point.
@@ -348,3 +489,69 @@ def _calculate_sde_parameters(beta0, beta1, beta3, beta4, v_mp, i_mp, v_oc):
348489
else: # I0_voc > 0
349490
I0 = I0_voc
350491
return (IL, I0, Rsh, Rs, nNsVth)
492+
493+
494+
def _system_of_equations_desoto(params, specs):
495+
"""Evaluates the systems of equations used to solve for the single
496+
diode equation parameters. Function designed to be used by
497+
scipy.optimize.root() in fit_sdm_desoto().
498+
499+
Parameters
500+
----------
501+
params: ndarray
502+
Array with parameters of the De Soto single diode model. Must be
503+
given in the following order: IL, Io, a, Rsh, Rs
504+
specs: tuple
505+
Specifications of pv module given by manufacturer. Must be given
506+
in the following order: Isc, Voc, Imp, Vmp, beta_oc, alpha_sc
507+
508+
Returns
509+
-------
510+
system of equations to solve with scipy.optimize.root().
511+
512+
513+
References
514+
----------
515+
[1] W. De Soto et al., "Improvement and validation of a model for
516+
photovoltaic array performance", Solar Energy, vol 80, pp. 78-88,
517+
2006.
518+
519+
[2] John A Duffie, William A Beckman, "Solar Engineering of Thermal
520+
Processes", Wiley, 2013
521+
"""
522+
523+
# six input known variables
524+
Isc, Voc, Imp, Vmp, beta_oc, alpha_sc, EgRef, dEgdT, Tref, k = specs
525+
526+
# five parameters vector to find
527+
IL, Io, a, Rsh, Rs = params
528+
529+
# five equation vector
530+
y = [0, 0, 0, 0, 0]
531+
532+
# 1st equation - short-circuit - eq(3) in [1]
533+
y[0] = Isc - IL + Io * np.expm1(Isc * Rs / a) + Isc * Rs / Rsh
534+
535+
# 2nd equation - open-circuit Tref - eq(4) in [1]
536+
y[1] = -IL + Io * np.expm1(Voc / a) + Voc / Rsh
537+
538+
# 3rd equation - Imp & Vmp - eq(5) in [1]
539+
y[2] = Imp - IL + Io * np.expm1((Vmp + Imp * Rs) / a) \
540+
+ (Vmp + Imp * Rs) / Rsh
541+
542+
# 4th equation - Pmp derivated=0 - eq23.2.6 in [2]
543+
# caution: eq(6) in [1] has a sign error
544+
y[3] = Imp \
545+
- Vmp * ((Io / a) * np.exp((Vmp + Imp * Rs) / a) + 1.0 / Rsh) \
546+
/ (1.0 + (Io * Rs / a) * np.exp((Vmp + Imp * Rs) / a) + Rs / Rsh)
547+
548+
# 5th equation - open-circuit T2 - eq (4) at temperature T2 in [1]
549+
T2 = Tref + 2
550+
Voc2 = (T2 - Tref) * beta_oc + Voc # eq (7) in [1]
551+
a2 = a * T2 / Tref # eq (8) in [1]
552+
IL2 = IL + alpha_sc * (T2 - Tref) # eq (11) in [1]
553+
Eg2 = EgRef * (1 + dEgdT * (T2 - Tref)) # eq (10) in [1]
554+
Io2 = Io * (T2 / Tref)**3 * np.exp(1 / k * (EgRef/Tref - Eg2/T2)) # eq (9)
555+
y[4] = -IL2 + Io2 * np.expm1(Voc2 / a2) + Voc2 / Rsh # eq (4) at T2
556+
557+
return y

pvlib/test/test_ivtools.py

+29
Original file line numberDiff line numberDiff line change
@@ -102,6 +102,35 @@ def test_fit_sdm_cec_sam(get_cec_params_cansol_cs5p_220p):
102102
cells_in_series=1, temp_ref=25)
103103

104104

105+
@requires_scipy
106+
def test_fit_sdm_desoto():
107+
result, _ = ivtools.fit_sdm_desoto(v_mp=31.0, i_mp=8.71, v_oc=38.3,
108+
i_sc=9.43, alpha_sc=0.005658,
109+
beta_voc=-0.13788,
110+
cells_in_series=60)
111+
result_expected = {'I_L_ref': 9.45232,
112+
'I_o_ref': 3.22460e-10,
113+
'a_ref': 1.59128,
114+
'R_sh_ref': 125.798,
115+
'R_s': 0.297814,
116+
'alpha_sc': 0.005658,
117+
'EgRef': 1.121,
118+
'dEgdT': -0.0002677,
119+
'irrad_ref': 1000,
120+
'temp_ref': 25}
121+
assert np.allclose(pd.Series(result), pd.Series(result_expected),
122+
rtol=1e-4)
123+
124+
125+
@requires_scipy
126+
def test_fit_sdm_desoto_failure():
127+
with pytest.raises(RuntimeError) as exc:
128+
ivtools.fit_sdm_desoto(v_mp=31.0, i_mp=8.71, v_oc=38.3, i_sc=9.43,
129+
alpha_sc=0.005658, beta_voc=-0.13788,
130+
cells_in_series=10)
131+
assert ('Parameter estimation failed') in str(exc.value)
132+
133+
105134
@pytest.fixture
106135
def get_bad_iv_curves():
107136
# v1, i1 produces a bad value for I0_voc

0 commit comments

Comments
 (0)