Skip to content

Commit c67c9ec

Browse files
authored
Add partial shading example using new reverse bias functionality (#968)
* Create plot_partial_module_shading_simple.py * add commentary * rename parameters to cell_parameters, stickler * remame combine to combine_series * add example cell curves and commentary * add module curve plot * lint * lint again, don't know how to fix this * updates from review * all tcell=25, remove method parameter * stickler
1 parent f6f9207 commit c67c9ec

File tree

1 file changed

+344
-0
lines changed

1 file changed

+344
-0
lines changed
Lines changed: 344 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,344 @@
1+
"""
2+
Calculating power loss from partial module shading
3+
==================================================
4+
5+
Example of modeling cell-to-cell mismatch loss from partial module shading.
6+
"""
7+
8+
# %%
9+
# Even though the PV cell is the primary power generation unit, PV modeling is
10+
# often done at the module level for simplicity because module-level parameters
11+
# are much more available and it significantly reduces the computational scope
12+
# of the simulation. However, module-level simulations are too coarse to be
13+
# able to model effects like cell to cell mismatch or partial shading. This
14+
# example calculates cell-level IV curves and combines them to reconstruct
15+
# the module-level IV curve. It uses this approach to find the maximum power
16+
# under various shading and irradiance conditions.
17+
#
18+
# The primary functions used here are:
19+
#
20+
# - :py:meth:`pvlib.pvsystem.calcparams_desoto` to estimate the single
21+
# diode equation parameters at some specified operating conditions.
22+
# - :py:meth:`pvlib.singlediode.bishop88` to calculate the full cell IV curve,
23+
# including the reverse bias region.
24+
#
25+
# .. note::
26+
#
27+
# This example requires the reverse bias functionality added in pvlib 0.7.2
28+
#
29+
# .. warning::
30+
#
31+
# Modeling partial module shading is complicated and depends significantly
32+
# on the module's electrical topology. This example makes some simplifying
33+
# assumptions that are not generally applicable. For instance, it assumes
34+
# that shading only applies to beam irradiance (*i.e.* all cells receive
35+
# the same amount of diffuse irradiance) and cell temperature is uniform
36+
# and not affected by cell-level irradiance variation.
37+
38+
from pvlib import pvsystem, singlediode
39+
import pandas as pd
40+
import numpy as np
41+
from scipy.interpolate import interp1d
42+
import matplotlib.pyplot as plt
43+
44+
from scipy.constants import e as qe, k as kB
45+
46+
# For simplicity, use cell temperature of 25C for all calculations.
47+
# kB is J/K, qe is C=J/V
48+
# kB * T / qe -> V
49+
Vth = kB * (273.15+25) / qe
50+
51+
cell_parameters = {
52+
'I_L_ref': 8.24,
53+
'I_o_ref': 2.36e-9,
54+
'a_ref': 1.3*Vth,
55+
'R_sh_ref': 1000,
56+
'R_s': 0.00181,
57+
'alpha_sc': 0.0042,
58+
'breakdown_factor': 2e-3,
59+
'breakdown_exp': 3,
60+
'breakdown_voltage': -15,
61+
}
62+
63+
# %%
64+
# Simulating a cell IV curve
65+
# --------------------------
66+
#
67+
# First, calculate IV curves for individual cells. The process is as follows:
68+
#
69+
# 1) Given a set of cell parameters at reference conditions and the operating
70+
# conditions of interest (irradiance and temperature), use a single-diode
71+
# model to calculate the single diode equation parameters for the cell at
72+
# the operating conditions. Here we use the De Soto model via
73+
# :py:func:`pvlib.pvsystem.calcparams_desoto`.
74+
# 2) The single diode equation cannot be solved analytically, so pvlib has
75+
# implemented a couple methods of solving it for us. However, currently
76+
# only the Bishop '88 method (:py:func:`pvlib.singlediode.bishop88`) has
77+
# the ability to model the reverse bias characteristic in addition to the
78+
# forward characteristic. Depending on the nature of the shadow, it is
79+
# sometimes necessary to model the reverse bias portion of the IV curve,
80+
# so we use the Bishop '88 method here. This gives us a set of (V, I)
81+
# points on the cell's IV curve.
82+
83+
84+
def simulate_full_curve(parameters, Geff, Tcell, ivcurve_pnts=1000):
85+
"""
86+
Use De Soto and Bishop to simulate a full IV curve with both
87+
forward and reverse bias regions.
88+
"""
89+
# adjust the reference parameters according to the operating
90+
# conditions using the De Soto model:
91+
sde_args = pvsystem.calcparams_desoto(
92+
Geff,
93+
Tcell,
94+
alpha_sc=parameters['alpha_sc'],
95+
a_ref=parameters['a_ref'],
96+
I_L_ref=parameters['I_L_ref'],
97+
I_o_ref=parameters['I_o_ref'],
98+
R_sh_ref=parameters['R_sh_ref'],
99+
R_s=parameters['R_s'],
100+
)
101+
# sde_args has values:
102+
# (photocurrent, saturation_current, resistance_series,
103+
# resistance_shunt, nNsVth)
104+
105+
# Use Bishop's method to calculate points on the IV curve with V ranging
106+
# from the reverse breakdown voltage to open circuit
107+
kwargs = {
108+
'breakdown_factor': parameters['breakdown_factor'],
109+
'breakdown_exp': parameters['breakdown_exp'],
110+
'breakdown_voltage': parameters['breakdown_voltage'],
111+
}
112+
v_oc = singlediode.bishop88_v_from_i(
113+
0.0, *sde_args, **kwargs
114+
)
115+
# ideally would use some intelligent log-spacing to concentrate points
116+
# around the forward- and reverse-bias knees, but this is good enough:
117+
vd = np.linspace(0.99*kwargs['breakdown_voltage'], v_oc, ivcurve_pnts)
118+
119+
ivcurve_i, ivcurve_v, _ = singlediode.bishop88(vd, *sde_args, **kwargs)
120+
return pd.DataFrame({
121+
'i': ivcurve_i,
122+
'v': ivcurve_v,
123+
})
124+
125+
126+
# %%
127+
# Now that we can calculate cell-level IV curves, let's compare a
128+
# fully-illuminated cell's curve to a shaded cell's curve. Note that shading
129+
# typically does not reduce a cell's illumination to zero -- tree shading and
130+
# row-to-row shading block the beam portion of irradiance but leave the diffuse
131+
# portion largely intact. In this example plot, we choose :math:`200 W/m^2`
132+
# as the amount of irradiance received by a shaded cell.
133+
134+
def plot_curves(dfs, labels, title):
135+
"""plot the forward- and reverse-bias portions of an IV curve"""
136+
fig, axes = plt.subplots(1, 2, sharey=True, figsize=(5, 3))
137+
for df, label in zip(dfs, labels):
138+
df.plot('v', 'i', label=label, ax=axes[0])
139+
df.plot('v', 'i', label=label, ax=axes[1])
140+
axes[0].set_xlim(right=0)
141+
axes[0].set_ylim([0, 25])
142+
axes[1].set_xlim([0, df['v'].max()*1.5])
143+
axes[0].set_ylabel('current [A]')
144+
axes[0].set_xlabel('voltage [V]')
145+
axes[1].set_xlabel('voltage [V]')
146+
fig.suptitle(title)
147+
fig.tight_layout()
148+
return axes
149+
150+
151+
cell_curve_full_sun = simulate_full_curve(cell_parameters, Geff=1000, Tcell=25)
152+
cell_curve_shaded = simulate_full_curve(cell_parameters, Geff=200, Tcell=25)
153+
ax = plot_curves([cell_curve_full_sun, cell_curve_shaded],
154+
labels=['Full Sun', 'Shaded'],
155+
title='Cell-level reverse- and forward-biased IV curves')
156+
157+
# %%
158+
# This figure shows how a cell's current decreases roughly in proportion to
159+
# the irradiance reduction from shading, but voltage changes much less.
160+
# At the cell level, the effect of shading is essentially to shift the I-V
161+
# curve down to lower currents rather than change the curve's shape.
162+
#
163+
# Note that the forward and reverse curves are plotted separately to
164+
# accommodate the different voltage scales involved -- a normal crystalline
165+
# silicon cell reaches only ~0.6V in forward bias, but can get to -10 to -20V
166+
# in reverse bias.
167+
#
168+
# Combining cell IV curves to create a module IV curve
169+
# ----------------------------------------------------
170+
#
171+
# To combine the individual cell IV curves and form a module's IV curve,
172+
# the cells in each substring must be added in series. The substrings are
173+
# in series as well, but with parallel bypass diodes to protect from reverse
174+
# bias voltages. To add in series, the voltages for a given current are
175+
# added. However, because each cell's curve is discretized and the currents
176+
# might not line up, we align each curve to a common set of current values
177+
# with interpolation.
178+
179+
180+
def interpolate(df, i):
181+
"""convenience wrapper around scipy.interpolate.interp1d"""
182+
f_interp = interp1d(np.flipud(df['i']), np.flipud(df['v']), kind='linear',
183+
fill_value='extrapolate')
184+
return f_interp(i)
185+
186+
187+
def combine_series(dfs):
188+
"""
189+
Combine IV curves in series by aligning currents and summing voltages.
190+
The current range is based on the first curve's current range.
191+
"""
192+
df1 = dfs[0]
193+
imin = df1['i'].min()
194+
imax = df1['i'].max()
195+
i = np.linspace(imin, imax, 1000)
196+
v = 0
197+
for df2 in dfs:
198+
v_cell = interpolate(df2, i)
199+
v += v_cell
200+
return pd.DataFrame({'i': i, 'v': v})
201+
202+
203+
# %%
204+
# Rather than simulate all 72 cells in the module, we'll assume that there
205+
# are only three types of cells (fully illuminated, fully shaded, and
206+
# partially shaded), and within each type all cells behave identically. This
207+
# means that simulating one cell from each type (for three cell simulations
208+
# total) is sufficient to model the module as a whole.
209+
#
210+
# This function also models the effect of bypass diodes in parallel with each
211+
# substring. Bypass diodes are normally inactive but conduct when substring
212+
# voltage becomes sufficiently negative, presumably due to the substring
213+
# entering reverse bias from mismatch between substrings. In that case the
214+
# substring's voltage is clamped to the diode's trigger voltage (assumed to
215+
# be 0.5V here).
216+
217+
def simulate_module(cell_parameters, poa_direct, poa_diffuse, Tcell,
218+
shaded_fraction, cells_per_string=24, strings=3):
219+
"""
220+
Simulate the IV curve for a partially shaded module.
221+
The shade is assumed to be coming up from the bottom of the module when in
222+
portrait orientation, so it affects all substrings equally.
223+
For simplicity, cell temperature is assumed to be uniform across the
224+
module, regardless of variation in cell-level illumination.
225+
Substrings are assumed to be "down and back", so the number of cells per
226+
string is divided between two columns of cells.
227+
"""
228+
# find the number of cells per column that are in full shadow
229+
nrow = cells_per_string // 2
230+
nrow_full_shade = int(shaded_fraction * nrow)
231+
# find the fraction of shade in the border row
232+
partial_shade_fraction = 1 - (shaded_fraction * nrow - nrow_full_shade)
233+
234+
df_lit = simulate_full_curve(
235+
cell_parameters,
236+
poa_diffuse + poa_direct,
237+
Tcell)
238+
df_partial = simulate_full_curve(
239+
cell_parameters,
240+
poa_diffuse + partial_shade_fraction * poa_direct,
241+
Tcell)
242+
df_shaded = simulate_full_curve(
243+
cell_parameters,
244+
poa_diffuse,
245+
Tcell)
246+
# build a list of IV curves for a single column of cells (half a substring)
247+
include_partial_cell = (shaded_fraction < 1)
248+
half_substring_curves = (
249+
[df_lit] * (nrow - nrow_full_shade - 1)
250+
+ ([df_partial] if include_partial_cell else []) # noqa: W503
251+
+ [df_shaded] * nrow_full_shade # noqa: W503
252+
)
253+
substring_curve = combine_series(half_substring_curves)
254+
substring_curve['v'] *= 2 # turn half strings into whole strings
255+
# bypass diode:
256+
substring_curve['v'] = substring_curve['v'].clip(lower=-0.5)
257+
# no need to interpolate since we're just scaling voltage directly:
258+
substring_curve['v'] *= strings
259+
return substring_curve
260+
261+
# %%
262+
# Now let's see how shade affects the IV curves at the module level. For this
263+
# example, the bottom 10% of the module is shaded. Assuming 12 cells per
264+
# column, that means one row of cells is fully shaded and another row is
265+
# partially shaded. Even though only 10% of the module is shaded, the
266+
# maximum power is decreased by roughly 80%!
267+
#
268+
# Note the effect of the bypass diodes. Without bypass diodes, operating the
269+
# shaded module at the same current as the fully illuminated module would
270+
# create a reverse-bias voltage of several hundred volts! However, the diodes
271+
# prevent the reverse voltage from exceeding 1.5V (three diodes at 0.5V each).
272+
273+
274+
kwargs = {
275+
'cell_parameters': cell_parameters,
276+
'poa_direct': 800,
277+
'poa_diffuse': 200,
278+
'Tcell': 25
279+
}
280+
module_curve_full_sun = simulate_module(shaded_fraction=0, **kwargs)
281+
module_curve_shaded = simulate_module(shaded_fraction=0.1, **kwargs)
282+
ax = plot_curves([module_curve_full_sun, module_curve_shaded],
283+
labels=['Full Sun', 'Shaded'],
284+
title='Module-level reverse- and forward-biased IV curves')
285+
286+
# %%
287+
# Calculating shading loss across shading scenarios
288+
# -------------------------------------------------
289+
#
290+
# Clearly the module-level IV-curve is strongly affected by partial shading.
291+
# This heatmap shows the module maximum power under a range of partial shade
292+
# conditions, where "diffuse fraction" refers to the ratio
293+
# :math:`poa_{diffuse} / poa_{global}` and "shaded fraction" refers to the
294+
# fraction of the module that receives only diffuse irradiance.
295+
296+
297+
def find_pmp(df):
298+
"""simple function to find Pmp on an IV curve"""
299+
return df.product(axis=1).max()
300+
301+
302+
# find Pmp under different shading conditions
303+
data = []
304+
for diffuse_fraction in np.linspace(0, 1, 11):
305+
for shaded_fraction in np.linspace(0, 1, 51):
306+
307+
df = simulate_module(cell_parameters,
308+
poa_direct=(1-diffuse_fraction)*1000,
309+
poa_diffuse=diffuse_fraction*1000,
310+
Tcell=25,
311+
shaded_fraction=shaded_fraction)
312+
data.append({
313+
'fd': diffuse_fraction,
314+
'fs': shaded_fraction,
315+
'pmp': find_pmp(df)
316+
})
317+
318+
results = pd.DataFrame(data)
319+
results['pmp'] /= results['pmp'].max() # normalize power to 0-1
320+
results_pivot = results.pivot('fd', 'fs', 'pmp')
321+
plt.figure()
322+
plt.imshow(results_pivot, origin='lower', aspect='auto')
323+
plt.xlabel('shaded fraction')
324+
plt.ylabel('diffuse fraction')
325+
xlabels = ["{:0.02f}".format(fs) for fs in results_pivot.columns[::5]]
326+
ylabels = ["{:0.02f}".format(fd) for fd in results_pivot.index]
327+
plt.xticks(range(0, 5*len(xlabels), 5), xlabels)
328+
plt.yticks(range(0, len(ylabels)), ylabels)
329+
plt.title('Module P_mp across shading conditions')
330+
plt.colorbar()
331+
plt.show()
332+
# use this figure as the thumbnail:
333+
# sphinx_gallery_thumbnail_number = 3
334+
335+
# %%
336+
# The heatmap makes a few things evident:
337+
#
338+
# - When diffuse fraction is equal to 1, there is no beam irradiance to lose,
339+
# so shading has no effect on production.
340+
# - When shaded fraction is equal to 0, no irradiance is blocked, so module
341+
# output does not change with the diffuse fraction.
342+
# - Under sunny conditions (diffuse fraction < 0.5), module output is
343+
# significantly reduced after just the first cell is shaded
344+
# (1/12 = ~8% shaded fraction).

0 commit comments

Comments
 (0)