def test_process_junction_pn():
    from solcore.analytic_solar_cells.depletion_approximation import identify_layers, identify_parameters
    from solcore import material
    from solcore.structure import Layer, Junction
    from solcore.state import State

    Na = np.power(10, np.random.uniform(22, 25))
    Nd = np.power(10, np.random.uniform(23, 26))

    options = State()
    options.T = np.random.uniform(250, 350)

    Lp = np.power(10, np.random.uniform(-9, -6))  # Diffusion length
    Ln = np.power(10, np.random.uniform(-9, -6))  # Diffusion length

    GaAs_n = material("GaAs")(Nd=Nd, hole_diffusion_length=Ln)
    GaAs_p = material("GaAs")(Na=Na, electron_diffusion_length=Lp)

    p_width = np.random.uniform(500, 1000)*1e-9
    n_width = np.random.uniform(3000, 5000)*1e-9

    test_junc = Junction([Layer(p_width, GaAs_p, role="emitter"), Layer(n_width, GaAs_n, role="base")])

    id_top, id_bottom, pRegion, nRegion, iRegion, pn_or_np = identify_layers(test_junc)
    xn, xp, xi, sn, sp, ln, lp, dn, dp, Nd_c, Na_c, ni, es = identify_parameters(test_junc, options.T, pRegion, nRegion, iRegion)

    ni_expect = GaAs_n.ni

    assert [id_top, id_bottom] == approx([0, 1])
    assert pn_or_np == 'pn'
    assert [xn, xp, xi, sn, sp, ln, lp, dn, dp, Nd_c, Na_c, ni, es] == approx([n_width, p_width, 0, 0, 0, Lp, Ln, GaAs_n.hole_mobility * kb * options.T / q,
                                    GaAs_p.electron_mobility * kb * options.T / q, Nd, Na, ni_expect, GaAs_n.permittivity])
def test_process_junction_set_in_junction():
    from solcore.analytic_solar_cells.depletion_approximation import identify_layers, identify_parameters
    from solcore import material
    from solcore.structure import Layer, Junction
    from solcore.state import State

    options = State()
    options.T = np.random.uniform(250, 350)

    Lp = np.power(10, np.random.uniform(-9, -6))  # Diffusion length
    Ln = np.power(10, np.random.uniform(-9, -6))  # Diffusion length

    sn = np.power(10, np.random.uniform(0, 3))
    sp = np.power(10, np.random.uniform(0, 3))

    se = np.random.uniform(1, 20) * vacuum_permittivity

    GaAs_n = material("GaAs")()
    GaAs_p = material("GaAs")()
    GaAs_i = material("GaAs")()

    p_width = np.random.uniform(500, 1000)
    n_width = np.random.uniform(3000, 5000)
    i_width = np.random.uniform(100, 300) * 1e-9

    mun = np.power(10, np.random.uniform(-5, 0))
    mup = np.power(10, np.random.uniform(-5, 0))

    Vbi = np.random.uniform(0, 3)

    test_junc = Junction([
        Layer(p_width, GaAs_p, role="emitter"),
        Layer(i_width, GaAs_i, role="intrinsic"),
        Layer(n_width, GaAs_n, role="base")
    ],
                         sn=sn,
                         sp=sp,
                         permittivity=se,
                         ln=Ln,
                         lp=Lp,
                         mup=mup,
                         mun=mun,
                         Vbi=Vbi)

    id_top, id_bottom, pRegion, nRegion, iRegion, pn_or_np = identify_layers(
        test_junc)
    xn, xp, xi, sn_c, sp_c, ln, lp, dn, dp, Nd_c, Na_c, ni, es = identify_parameters(
        test_junc, options.T, pRegion, nRegion, iRegion)

    ni_expect = GaAs_n.ni

    assert [id_top, id_bottom] == approx([0, 2])
    assert pn_or_np == 'pn'
    assert [xn, xp, xi, sn_c, sp_c, ln, lp, dn, dp, Nd_c, Na_c, ni,
            es] == approx([
                n_width, p_width, i_width, sn, sp, Ln, Lp,
                mun * kb * options.T / q, mup * kb * options.T / q, 1, 1,
                ni_expect, se
            ])
Пример #3
0
def merge_dicts(*dict_args):
    """
    Given any number of dicts, shallow copy and merge into a new dict,
    precedence goes to key value pairs in latter dicts.
    """
    result = State()
    for dictionary in dict_args:
        result.update(dictionary)
    return result
def test_BL_correction():

    wl = np.linspace(800, 950, 4) * 1e-9

    GaAs = material("GaAs")()

    thick_cell = SolarCell([Layer(material=GaAs, width=si("20um"))])

    opts = State()
    opts.position = None
    prepare_solar_cell(thick_cell, opts)
    position = np.arange(0, thick_cell.width, 1e-9)
    opts.position = position
    opts.recalculate_absorption = True
    opts.no_back_reflexion = False

    opts.BL_correction = False
    opts.wavelength = wl
    solve_tmm(thick_cell, opts)

    no_corr = thick_cell.absorbed

    opts.BL_correction = True

    solve_tmm(thick_cell, opts)

    with_corr = thick_cell.absorbed

    assert with_corr == approx(
        np.array([6.71457872e-01, 6.75496354e-01, 2.09738887e-01, 0]))
    assert no_corr == approx(
        np.array([6.71457872e-01, 6.75496071e-01, 2.82306407e-01, 0]))
Пример #5
0
def qe_pdd(junction, options):
    """ Calculates the quantum efficiency of the device at short circuit. Internally it calls ShortCircuit

    :param device: A dictionary containing the device structure. See PDD.DeviceStructure
    :param wavelengths: Array with the wavelengths (in m)
    :param photon_flux: Array with the photon_flux (in photons/m2/m) corresponding to the above wavelengths
    :param output_info: Indicates how much information must be printed by the fortran solver (0=min, 2=max)
    :param rs: Series resistance. Default=0
    :return: The internal and external quantum efficiencies, in adition to the output of ShortCircuit.
    """
    print('Solving quantum efficiency...')
    output_info = options.output_qe

    short_circuit_pdd(junction, options)

    dd.runiqe(output_info)
    print('...done!\n')

    output = {}
    output['QE'] = DumpQE()
    output['QE']['wavelengths'] = options.wavelength

    junction.qe_data = State(**output)

    # The EQE is actually the IQE inside the fortran solver due to an error in the naming --> to be changed
    junction.eqe = interp1d(options.wavelength,
                            output['QE']['IQE'],
                            kind='linear',
                            bounds_error=False,
                            assume_sorted=True,
                            fill_value=(output['QE']['IQE'][0],
                                        output['QE']['IQE'][-1]))
Пример #6
0
def equilibrium_pdd(junction, options):
    """ Solves the PDD equations under equilibrium: in the dark with no external current and zero applied voltage.

    :param junction: A junction object
    :param options: Options to be passed to the solver
    :return: None
    """
    T = options.T
    wl = options.wavelength
    output_info = options.output_equilibrium

    SetConvergenceParameters(options)
    SetMeshParameters(options)
    SetRecombinationParameters(options)

    device = CreateDeviceStructure('Junction', T=T, layers=junction)

    print('Solving equilibrium...')
    output = ProcessStructure(device, options.meshpoints, wavelengths=wl)
    dd.gen = 0

    dd.equilibrium(output_info)
    print('...done!\n')

    output['Bandstructure'] = DumpBandStructure()
    junction.equilibrium_data = State(**output)
Пример #7
0
def qe_pdd(junction, options):
    """ Calculates the quantum efficiency of the device at short circuit. Internally it calls ShortCircuit

    :param junction: A junction object
    :param options: Options to be passed to the solver
    :return: None
    """
    print('Solving quantum efficiency...')
    output_info = options.output_qe

    short_circuit_pdd(junction, options)

    dd.runiqe(output_info)
    print('...done!\n')

    output = {}
    output['QE'] = DumpQE()
    output['QE']['WL'] = options.wavelength

    z = junction.short_circuit_data['Bandstructure']['x']
    absorbed_per_wl = np.trapz(junction.absorbed(z), z, axis=0)

    # This is redundant but useful to keep the same format than the other solvers
    junction.qe = State(**output['QE'])
    junction.qe['IQE'] = junction.qe['EQE'] / np.maximum(
        absorbed_per_wl, 0.00001)

    # The EQE is actually the IQE inside the fortran solver due to an error in the naming --> to be changed
    junction.eqe = interp1d(options.wavelength,
                            output['QE']['EQE'],
                            kind='linear',
                            bounds_error=False,
                            assume_sorted=True,
                            fill_value=(output['QE']['EQE'][0],
                                        output['QE']['EQE'][-1]))
Пример #8
0
def short_circuit_pdd(junction, options):
    """ Solves the devices electronic properties at short circuit. Internally, it calls Equilibrium.

    :param junction: A junction object
    :param options: Options to be passed to the solver
    :return: None
    """

    # We run equilibrium
    equilibrium_pdd(junction, options)

    wl = options.wavelength
    wl, ph = options.light_source.spectrum(x=wl,
                                           output_units='photon_flux_per_m')

    z = junction.equilibrium_data['Bandstructure']['x']
    absorption = junction.absorbed if hasattr(junction,
                                              'absorbed') else lambda x: 0
    abs_profile = CalculateAbsorptionProfile(z, wl, absorption)

    dd.set_generation(abs_profile)
    dd.gen = 1

    dd.illumination(ph)

    dd.lightsc(options.output_sc, 1)
    print('...done!\n')

    output = {'Bandstructure': DumpBandStructure()}
    junction.short_circuit_data = State(**output)
def qe_detailed_balance(junction, wl):
    """ Calculates the EQE and the IQE in the detailed balanced limit... which is equal to the absorptance and the absorbed fraction of the light since, by definition of detailed balanced, the carrier collection is perfect.

    :param junction: Junction object of kind "DB"
    :param wl: wavelength in m
    :return: None
    """
    absorptance_detailed_balance(junction)

    # This is true just in the detailed balance limit
    junction.iqe = junction.absorptance

    try:
        z = np.linspace(0, junction.width, 1001)
        all_abs = junction.absorbed(z)
        abs_vs_wl = np.trapz(all_abs, z, axis=0)
        junction.eqe = interp1d(wl,
                                abs_vs_wl,
                                bounds_error=False,
                                fill_value=(0, 0))

    except AttributeError:
        junction.eqe = junction.absorptance

    junction.qe = State({
        'WL': wl,
        'IQE': junction.iqe(wl),
        'EQE': junction.eqe(wl)
    })
Пример #10
0
def qe_pdd(junction, options):
    """ Calculates the quantum efficiency of the device at short circuit. Internally it calls ShortCircuit

    :param junction: A junction object
    :param options: Options to be passed to the solver
    :return: None
    """
    print('Solving quantum efficiency...')
    output_info = options.output_qe

    short_circuit_pdd(junction, options)

    dd.runiqe(output_info)
    print('...done!\n')

    output = {}
    output['QE'] = DumpQE()
    output['QE']['wavelengths'] = options.wavelength

    junction.qe_data = State(**output['QE'])

    # The EQE is actually the IQE inside the fortran solver due to an error in the naming --> to be changed
    junction.eqe = interp1d(options.wavelength,
                            output['QE']['EQE'],
                            kind='linear',
                            bounds_error=False,
                            assume_sorted=True,
                            fill_value=(output['QE']['EQE'][0],
                                        output['QE']['EQE'][-1]))
Пример #11
0
def short_circuit_pdd(junction, options):
    """ Solves the devices electronic properties at short circuit. Internally, it calls Equilibrium.

    :param device: A dictionary containing the device structure. See PDD.DeviceStructure
    :param output_info: Indicates how much information must be printed by the fortran solver (0=min, 2=max)
    :param rs: Series resistance. Default=0
    :param wavelengths: Array with the wavelengths (in m)
    :param photon_flux: Array with the photon_flux (in photons/m2/m) corresponding to the above wavelengths
    :return: A dictionary containing the device properties as a function of the position at short circuit.
    """

    # We run equilibrium
    equilibrium_pdd(junction, options)

    wl = options.wavelength
    wl, ph = options.light_source.spectrum(x=wl,
                                           output_units='photon_flux_per_m')

    z = junction.equilibrium_data['Bandstructure']['x']
    absorption = junction.absorbed if hasattr(junction,
                                              'absorbed') else lambda x: 0
    abs_profile = CalculateAbsorptionProfile(z, wl, absorption)

    dd.set_generation(abs_profile)
    dd.gen = 1

    dd.illumination(ph)

    dd.lightsc(options.output_sc, 1)
    print('...done!\n')

    output = {'Bandstructure': DumpBandStructure()}
    junction.short_circuit_data = State(**output)
Пример #12
0
def equilibrium_pdd(junction, options):
    """ Solves the Poisson-DD equations under equilibrium: in the dark with no external current and zero applied voltage. Internally, it calls *ProcessStructure*. Absorption coeficients are not calculated unless *wavelengths* is given as input.

    :param device: A dictionary containing the device structure. See PDD.DeviceStructure
    :param output_info: Indicates how much information must be printed by the fortran solver (1=less, 2=more)
    :param wavelengths: (Optional) Wavelengths at which to calculate the optical properties.
    :return: Dictionary containing the device properties as a function of the position at equilibrium.
    """
    T = options.T
    wl = options.wavelength
    output_info = options.output_equilibrium

    SetConvergenceParameters(options)
    SetMeshParameters(options)
    SetRecombinationParameters(options)

    device = CreateDeviceStructure('Junction', T=T, layers=junction)

    print('Solving equilibrium...')
    output = ProcessStructure(device, options.meshpoints, wavelengths=wl)
    dd.gen = 0

    dd.equilibrium(output_info)
    print('...done!\n')

    output['Bandstructure'] = DumpBandStructure()
    junction.equilibrium_data = State(**output)
Пример #13
0
def da_options():
    from solcore.state import State

    options = State()
    wl = np.linspace(290, 700, 150) * 1e-9
    options.T = np.random.uniform(250, 350)
    options.wavelength = wl
    options.light_source = da_light_source()
    options.position = None
    options.internal_voltages = np.linspace(-6, 4, 20)
    options.da_mode = 'bvp'

    return options
Пример #14
0
def test_inc_coh_tmm():
    GaInP = material("GaInP")(In=0.5)
    GaAs = material("GaAs")()
    Ge = material("Ge")()

    optical_struct = SolarCell([
        Layer(material=GaInP, width=si("5000nm")),
        Layer(material=GaAs, width=si("200nm")),
        Layer(material=GaAs, width=si("5um")),
        Layer(material=Ge, width=si("50um")),
    ])

    wl = np.linspace(400, 1200, 5) * 1e-9

    options = State()
    options.wavelength = wl
    options.optics_method = "TMM"
    options.no_back_reflection = False
    options.BL_correction = True
    options.recalculate_absorption = True

    c_list = [
        ["c", "c", "c", "c"],
        ["c", "c", "c", "i"],
        ["c", "i", "i", "c"],
        ["i", "i", "i", "i"],
    ]

    results = []
    for i1, cl in enumerate(c_list):
        options.coherency_list = cl
        solar_cell_solver(optical_struct, "optics", options)
        results.append(optical_struct.absorbed)

    A_calc = np.stack(results)
    A_data = np.array(
        [[0.5742503, 0.67956899, 0.73481184, 0.725372, 0.76792856],
         [0.5742503, 0.67956899, 0.73481184, 0.725372, 0.76792856],
         [0.5742503, 0.67956899, 0.73474943, 0.70493469, 0.70361194],
         [0.5742503, 0.67956899, 0.70927724, 0.71509221, 0.71592772]])
    assert A_calc == approx(A_data)
Пример #15
0
def iv_pdd(junction, options):
    """ Calculates the IV curve of the device between 0 V and a given voltage. Depending on the options, the IV will be calculated in the dark (calling the equilibrium_pdd function) or under illumination (calling the short_circuit_pdd function). If the voltage range has possitive and negative values, the problem is solved twice: from 0 V to the maximu positive and from 0 V to the maximum negative, concatenating the results afterwards.

    :param junction: A junction object
    :param options: Options to be passed to the solver
    :return: None
    """
    print('Solving IV...')
    light = options.light_iv
    output_info = options.output_iv
    T = options.T

    junction.voltage = options.internal_voltages

    Eg = 10
    for layer in junction:
        Eg = min([Eg, layer.material.band_gap / q])

    s = True if junction[0].material.Na >= junction[0].material.Nd else False

    if s:
        vmax = min(Eg + 3 * kb * T / q, max(junction.voltage))
        vmin = min(junction.voltage)
    else:
        vmax = max(junction.voltage)
        vmin = max(-Eg - 3 * kb * T / q, min(junction.voltage))

    vstep = junction.voltage[1] - junction.voltage[0]

    # POSITIVE RANGE
    output_pos = None
    if vmax > 0:
        if light:
            short_circuit_pdd(junction, options)
        else:
            equilibrium_pdd(junction, options)

        dd.runiv(vmax, vstep, output_info, 0)
        output_pos = {'Bandstructure': DumpBandStructure(), 'IV': DumpIV()}

    # NEGATIVE RANGE
    output_neg = None
    if vmin < 0:
        if light:
            short_circuit_pdd(junction, options)
        else:
            equilibrium_pdd(junction, options)

        dd.runiv(vmin, (-1 * vstep), output_info, 0)
        output_neg = {'Bandstructure': DumpBandStructure(), 'IV': DumpIV()}

    print('...done!\n')

    # Now we need to put together the data for the possitive and negative regions.
    junction.pdd_data = State({
        'possitive_V': output_pos,
        'negative_V': output_neg
    })

    if output_pos is None:
        V = output_neg['IV']['V']
        J = output_neg['IV']['J']
        Jrad = output_neg['IV']['Jrad']
        Jsrh = output_neg['IV']['Jsrh']
        Jaug = output_neg['IV']['Jaug']
        Jsur = output_neg['IV']['Jsur']
    elif output_neg is None:
        V = output_pos['IV']['V']
        J = output_pos['IV']['J']
        Jrad = output_pos['IV']['Jrad']
        Jsrh = output_pos['IV']['Jsrh']
        Jaug = output_pos['IV']['Jaug']
        Jsur = output_pos['IV']['Jsur']
    else:
        V = np.concatenate(
            (output_neg['IV']['V'][:0:-1], output_pos['IV']['V']))
        J = np.concatenate(
            (output_neg['IV']['J'][:0:-1], output_pos['IV']['J']))
        Jrad = np.concatenate(
            (output_neg['IV']['Jrad'][:0:-1], output_pos['IV']['Jrad']))
        Jsrh = np.concatenate(
            (output_neg['IV']['Jsrh'][:0:-1], output_pos['IV']['Jsrh']))
        Jaug = np.concatenate(
            (output_neg['IV']['Jaug'][:0:-1], output_pos['IV']['Jaug']))
        Jsur = np.concatenate(
            (output_neg['IV']['Jsur'][:0:-1], output_pos['IV']['Jsur']))

    # Finally, we calculate the currents at the desired voltages
    R_shunt = min(junction.R_shunt, 1e14) if hasattr(junction,
                                                     'R_shunt') else 1e14

    junction.current = np.interp(junction.voltage, V,
                                 J) + junction.voltage / R_shunt
    Jrad = np.interp(junction.voltage, V, Jrad)
    Jsrh = np.interp(junction.voltage, V, Jsrh)
    Jaug = np.interp(junction.voltage, V, Jaug)
    Jsur = np.interp(junction.voltage, V, Jsur)

    junction.iv = interp1d(junction.voltage,
                           junction.current,
                           kind='linear',
                           bounds_error=False,
                           assume_sorted=True,
                           fill_value=(junction.current[0],
                                       junction.current[-1]))
    junction.recombination_currents = State({
        "Jrad": Jrad,
        "Jsrh": Jsrh,
        "Jaug": Jaug,
        "Jsur": Jsur
    })
import numpy as np
from scipy.interpolate import interp1d
from scipy.integrate import solve_bvp, quad_vec
from functools import partial

from solcore.constants import kb, q
from solcore.science_tracker import science_reference
from solcore.state import State
from solcore.light_source import LightSource

da_options = State()
da_options.da_mode = 'bvp'


def identify_layers(junction):
    # First we have to figure out if we are talking about a PN, NP, PIN or NIP junction
    # We search for the emitter and check if it is n-type or p-type
    idx = 0
    pn_or_np = 'pn'
    homojunction = True

    for layer in junction:
        if layer.role.lower() != 'emitter':
            idx += 1
        else:
            Na = 0
            Nd = 0
            if hasattr(layer.material, 'Na'): Na = layer.material.Na
            if hasattr(layer.material, 'Nd'): Nd = layer.material.Nd
            if Na < Nd:
                pn_or_np = "np"
Пример #17
0
from solcore import asUnit, constants
from solcore.state import State
from .ddModel import driftdiffusion as dd
from .DeviceStructure import LoadAbsorption, CreateDeviceStructure, CalculateAbsorptionProfile

Epsi0 = constants.vacuum_permittivity
q = constants.q
pi = constants.pi
h = constants.h
kb = constants.kb
m0 = constants.electron_mass

log = dd.log_file

pdd_options = State()

# Mesh control
pdd_options.meshpoints = -400
pdd_options.growth_rate = 0.7
pdd_options.coarse = 20e-9
pdd_options.fine = 1e-9
pdd_options.ultrafine = 0.2e-9

# Convergence control
pdd_options.clamp = 20
pdd_options.nitermax = 100
pdd_options.ATol = 1e-14
pdd_options.RTol = 1e-6

# Recombination control
Пример #18
0
def iv_depletion(junction, options):
    """ Calculates the IV curve of a junction object using the depletion approximation as described in J. Nelson, “The Physics of Solar Cells”, Imperial College Press (2003). The junction is then updated with an "iv" function that calculates the IV curve at any voltage.

    :param junction: A junction object.
    :param options: Solver options.
    :return: None.
    """

    science_reference(
        'Depletion approximation',
        'J. Nelson, “The Physics of Solar Cells”, Imperial College Press (2003).'
    )

    junction.voltage = options.internal_voltages
    T = options.T

    # First we have to figure out if we are talking about a PN, NP, PIN or NIP junction
    sn = 0 if not hasattr(junction, "sn") else junction.sn
    sp = 0 if not hasattr(junction, "sp") else junction.sp

    # We search for the emitter and check if it is n-type or p-type
    idx = 0
    pn_or_np = 'pn'

    for layer in junction:
        if layer.role is not 'emitter':
            idx += 1
        else:
            Na = 0
            Nd = 0
            if hasattr(layer.material, 'Na'): Na = layer.material.Na
            if hasattr(layer.material, 'Nd'): Nd = layer.material.Nd
            if Na < Nd:
                pn_or_np = "np"
                nRegion = junction[idx]
            else:
                pRegion = junction[idx]

            id_top = idx

            break

    # Now we check for an intrinsic region and, if there is, for the base.
    if junction[idx + 1].role is 'intrinsic':
        iRegion = junction[idx + 1]

        if junction[idx + 2].role is 'base':
            if pn_or_np == "pn":
                nRegion = junction[idx + 2]
                nidx = idx + 2
            else:
                pRegion = junction[idx + 2]
                pidx = idx + 2

            id_bottom = idx + 2

        else:
            raise RuntimeError(
                'ERROR processing junctions: A layer following the "intrinsic" layer must be defined as '
                '"base".')

    # If there is no intrinsic region, we check directly the base
    elif junction[idx + 1].role is 'base':
        if pn_or_np == "pn":
            nRegion = junction[idx + 1]
            nidx = idx + 1
        else:
            pRegion = junction[idx + 1]
            pidx = idx + 1
        iRegion = None

        id_bottom = idx + 1

    else:
        raise RuntimeError(
            'ERROR processing junctions: A layer following the "emitter" must be defined as "intrinsic"'
            'or "base".')

    # With all regions identified, it's time to start doing calculations
    kbT = kb * T
    Egap = nRegion.material.band_gap
    R_shunt = min(junction.R_shunt, 1e14) if hasattr(junction,
                                                     'R_shunt') else 1e14

    xp = pRegion.width
    xn = nRegion.width
    xi = 0 if iRegion is None else iRegion.width

    # Now we have to get all the material parameters needed for the calculation
    if hasattr(junction, "dielectric_constant"):
        es = junction.dielectric_constant
    else:
        es = nRegion.material.permittivity * vacuum_permittivity  # equal for n and p.  I hope.

    # For the diffusion lenght, subscript n and p refer to the carriers, electrons and holes
    if hasattr(junction, "ln"):
        ln = junction.ln
    else:
        ln = pRegion.material.electron_diffusion_length

    if hasattr(junction, "lp"):
        lp = junction.lp
    else:
        lp = nRegion.material.hole_diffusion_length

    # For the diffusion coefficient, n and p refer to the regions, n side and p side. Yeah, it's confusing...
    if hasattr(junction, "mup"):
        dp = junction.mup * kbT / q
    else:
        dp = pRegion.material.electron_mobility * kbT / q

    if hasattr(junction, "mun"):
        dn = junction.mun * kbT / q
    else:
        dn = nRegion.material.hole_mobility * kbT / q

    # Effective masses and effective density of states
    mEff_h = nRegion.material.eff_mass_hh_z * electron_mass
    mEff_e = pRegion.material.eff_mass_electron * electron_mass

    Nv = 2 * (mEff_h * kb * T / (2 * pi * hbar**2))**1.5  # Jenny p58
    Nc = 2 * (mEff_e * kb * T / (2 * pi * hbar**2))**1.5
    niSquared = Nc * Nv * np.exp(-Egap / (kb * T))
    ni = np.sqrt(niSquared)

    Na = pRegion.material.Na
    Nd = nRegion.material.Nd
    Vbi = (kbT / q) * np.log(Nd * Na / niSquared) if not hasattr(
        junction, "Vbi") else junction.Vbi  # Jenny p146

    # And now we account for the possible applied voltage, which can be, at most, equal to Vbi
    V = np.where(junction.voltage < Vbi - 0.001, junction.voltage, Vbi - 0.001)

    # It's time to calculate the depletion widths
    if not hasattr(junction, "wp") or not hasattr(junction, "wn"):

        if hasattr(
                junction, "depletion_approximation"
        ) and junction.depletion_approximation == "one-sided abrupt":
            print(
                "using one-sided abrupt junction approximation for depletion width"
            )
            science_reference(
                "Sze abrupt junction approximation",
                "Sze: The Physics of Semiconductor Devices, 2nd edition, John Wiley & Sons, Inc (2007)"
            )
            wp = np.sqrt(2 * es * (Vbi - V) / (q * Na))
            wn = np.sqrt(2 * es * (Vbi - V) / (q * Nd))

        else:
            wn = (-xi + np.sqrt(xi**2 + 2. * es * (Vbi - V) / q *
                                (1 / Na + 1 / Nd))) / (1 + Nd / Na)
            wp = (-xi + np.sqrt(xi**2 + 2. * es * (Vbi - V) / q *
                                (1 / Na + 1 / Nd))) / (1 + Na / Nd)

    wn = wn if not hasattr(junction, "wn") else junction.wn
    wp = wp if not hasattr(junction, "wp") else junction.wp

    w = wn + wp + xi

    # Now it is time to calculate currents
    if pn_or_np == "pn":
        l_top, l_bottom = ln, lp
        x_top, x_bottom = xp, xn
        w_top, w_bottom = wp, wn
        s_top, s_bottom = sp, sn
        d_top, d_bottom = dp, dn
        min_top, min_bot = niSquared / Na, niSquared / Nd
    else:
        l_bottom, l_top = ln, lp
        x_bottom, x_top = xp, xn
        w_bottom, w_top = wp, wn
        s_bottom, s_top = sp, sn
        d_bottom, d_top = dp, dn
        min_bot, min_top = niSquared / Na, niSquared / Nd

    JtopDark = get_j_top(x_top, w_top, l_top, s_top, d_top, V, min_top, T)
    JbotDark = get_j_bot(x_bottom, w_bottom, l_bottom, s_bottom, d_bottom, V,
                         min_bot, T)

    # hereby we define the subscripts to refer to the layer in which the current is generated:
    if pn_or_np == "pn":
        JnDark, JpDark = JbotDark, JtopDark
    else:
        JpDark, JnDark = JbotDark, JtopDark

    # These might not be the right lifetimes. Actually, they are not as they include all recombination processes, not
    # just SRH recombination, which is what the equation in Jenny, p159 refers to. Let´ leave them, for now.
    lifetime_n = ln**2 / dn
    lifetime_p = lp**2 / dp  # Jenny p163

    # Here we use the full version of the SRH recombination term as calculated by Sah et al. Works for positive bias
    # and moderately negative ones.
    science_reference(
        'SRH current term.',
        'C. T. Sah, R. N. Noyce, and W. Shockley, “Carrier Generation and Recombination in P-N Junctions and P-N Junction Characteristics,” presented at the Proceedings of the IRE, 1957, vol. 45, no. 9, pp. 1228–1243.'
    )
    Jrec = get_Jsrh(ni, V, Vbi, lifetime_p, lifetime_n, w, kbT)

    J_sc_top = 0
    J_sc_bot = 0
    J_sc_scr = 0

    if options.light_iv:

        widths = []
        for layer in junction:
            widths.append(layer.width)

        cum_widths = np.cumsum([0] + widths)

        g = junction.absorbed
        wl = options.wavelength
        wl_sp, ph = options.light_source.spectrum(
            output_units='photon_flux_per_m', x=wl)
        id_v0 = np.argmin(abs(V))

        # The contribution from the Emitter (top side).
        xa = cum_widths[id_top]
        xb = cum_widths[id_top + 1] - w_top[id_v0]
        deriv = get_J_sc_diffusion(xa,
                                   xb,
                                   g,
                                   d_top,
                                   l_top,
                                   min_top,
                                   s_top,
                                   wl,
                                   ph,
                                   side='top')
        J_sc_top = q * d_top * abs(deriv)

        # The contribution from the Base (bottom side).
        xa = cum_widths[id_bottom] + w_bottom[id_v0]
        xb = cum_widths[id_bottom + 1]
        deriv = get_J_sc_diffusion(xa,
                                   xb,
                                   g,
                                   d_bottom,
                                   l_bottom,
                                   min_bot,
                                   s_bottom,
                                   wl,
                                   ph,
                                   side='bottom')
        J_sc_bot = q * d_bottom * abs(deriv)

        # The contribution from the SCR (includes the intrinsic region, if present).
        xa = cum_widths[id_top + 1] - w_top[id_v0]
        xb = cum_widths[id_bottom] + w_bottom[id_v0]
        J_sc_scr = q * get_J_sc_SCR(xa, xb, g, wl, ph)

    # And, finally, we output the currents
    junction.current = Jrec + JnDark + JpDark + V / R_shunt - J_sc_top - J_sc_bot - J_sc_scr
    junction.iv = interp1d(junction.voltage,
                           junction.current,
                           kind='linear',
                           bounds_error=False,
                           assume_sorted=True,
                           fill_value=(junction.current[0],
                                       junction.current[-1]))
    junction.region_currents = State({
        "Jn_dif": JnDark,
        "Jp_dif": JpDark,
        "Jscr_srh": Jrec,
        "J_sc_top": J_sc_top,
        "J_sc_bot": J_sc_bot,
        "J_sc_scr": J_sc_scr
    })
Пример #19
0
from solcore.structure import Layer, Junction, TunnelJunction
from solcore.absorption_calculator import calculate_rat_rcwa, calculate_absorption_profile_rcwa

import numpy as np
import types
from solcore.state import State

rcwa_options = State()
rcwa_options.size = [500, 500]
rcwa_options.orders = 4
rcwa_options.theta = 0
rcwa_options.phi = 0
rcwa_options.pol = 'u'


def solve_rcwa(solar_cell, options):
    """ Calculates the reflection, transmission and absorption of a solar cell object using the transfer matrix method

    :param solar_cell:
    :param options:
    :return:
    """
    wl = options.wavelength

    # We include the shadowing losses
    initial = (1 - solar_cell.shading) if hasattr(solar_cell, 'shading') else 1

    # Now we calculate the absorbed and transmitted light. We first get all the relevant parameters from the objects
    widths = []
    offset = 0
    all_layers = []
Пример #20
0
def strain_calculation_parameters(substrate_material,
                                  layer_material,
                                  should_print=False,
                                  SO=False):
    """Will extract materials parameters and perform a bit of conditioning to the values.
    
    Returns a solcore State object with the following keys:
    
        -- ac, the conduction band deformation potential
        -- av, the valence band deformation potential
        -- b,  the valence band splitting deformation potential
        -- C11, element of the elastic stiffness tensor
        -- C12, element of the elastic stiffness tensor
        -- a0, the unstrained substrate lattice constant
        -- a, the unstrained layer lattice constant
        -- epsilon, in-plain strain
        -- epsilon_perp, out-plain strain
        -- e_xx, in-plain strain (e_xx = epsilon)
        -- e_yy, in-plain strain (e_yy = epsilon)
        -- e_zz, out-plain strain (e_zz = epsilon_perp)
        -- Tre, element of come matrix (Tre = e_xx + e_yy + e_zz)
        -- Pe, parameter use by Chuang
        -- Qe, parameter use by Chuang
        -- cb_hydrostatic_shift, CB moved by this amount
        -- vb_hydrostatic_shift, VB moved by this amount
        -- vb_shear_splitting, VB split by this amount (i.e. HH/LH separation)
        -- delta_Ec, final conduction band shift
        -- delta_Elh, final light hole band shift
        -- delta_Ehh, final heavy hole band shift
    
    Care has to be taken when calculating the energy shifts because different 
    sign conversion are used by different authors. Here we use the approach of
    S. L. Chuang, 'Physics of optoelectronic devices'. 
    
    Note that this requires that the 'a_v' deformation potential to be 
    positive, where as Vurgaftman defines this a negative!
    """

    sub = substrate_material
    mat = layer_material
    k = State()

    # deformation potentials
    k.av = abs(mat.a_v)  # make sure that av is positive for this calculation
    k.ac = mat.a_c  # Vurgaftman uses the convention that this is negative
    k.b = mat.b

    # Matrix elements from elastic stiffness tensor
    k.C11 = mat.c11
    k.C12 = mat.c12
    if should_print: print(sub, mat)
    # Strain fractions
    k.a0 = sub.lattice_constant
    k.a = mat.lattice_constant
    k.epsilon = (k.a0 - k.a) / k.a  # in-plain strain
    k.epsilon_perp = -2 * k.C12 / k.C11 * k.epsilon  # out-plain
    k.e_xx = k.epsilon
    k.e_yy = k.epsilon
    k.e_zz = k.epsilon_perp
    k.Tre = (k.e_xx + k.e_yy + k.e_zz)
    k.Pe = -k.av * k.Tre
    k.Qe = -k.b / 2 * (k.e_xx + k.e_yy - 2 * k.e_zz)

    k.cb_hydrostatic_shift = k.ac * k.Tre
    k.vb_hydrostatic_shift = k.av * k.Tre
    k.vb_shear_splitting = 2 * k.b * (1 + 2 * k.C12 / k.C11) * k.epsilon

    # Shifts and splittings
    k.delta_Ec = k.ac * k.Tre

    if should_print: print(k.ac, k.Tre)

    k.delta_Ehh = -k.Pe - k.Qe
    k.delta_Elh = -k.Pe + k.Qe
    k.delta_Eso = 0.0

    if SO:
        k.delta = mat.spin_orbit_splitting
        shift = k.delta**2 + 2 * k.delta * k.Qe + 9 * k.Qe**2
        k.delta_Elh = -k.Pe + 0.5 * (k.Qe - k.delta + np.sqrt(shift))
        k.delta_Eso = -k.Pe + 0.5 * (k.Qe - k.delta - np.sqrt(shift))

    strain_calculation_asserts(k, should_print=should_print)

    if should_print:
        print()
        print("Lattice:")
        print("a0", k.a0)
        print("a", k.a)
        print()
        print("Deformation potentials:")
        print("ac = ", solcore.asUnit(k.ac, 'eV'))
        print("av = ", solcore.asUnit(k.av, 'eV'))
        print("ac - av = ", solcore.asUnit(k.ac - k.av, 'eV'))
        print("b = ", solcore.asUnit(k.b, 'eV'))
        print()
        print("Matrix elements from elastic stiffness tensor:")
        print("C_11 = ", solcore.asUnit(k.C11, "GPa"))
        print("C_12 = ", solcore.asUnit(k.C12, "GPa"))
        print()
        print("Strain fractions:")
        print("e_xx = e_yy = epsilon = ", k.epsilon)
        print("e_zz = epsilon_perp = ", k.epsilon_perp)
        print("e_xx + e_yy + e_zz = Tre = ", k.Tre)
        print()
        print("Shifts and splittings:")
        print("Pe = -av * Tre = ", solcore.asUnit(k.Pe, 'eV'))
        print("Qe = -b/2*(e_xx + e_yy - 2*e_zz) = ",
              solcore.asUnit(k.Qe, 'eV'))
        print("dEc = ac * Tre = ", solcore.asUnit(k.delta_Ec, 'eV'))
        print("dEhh = av * Tre + b[1 + 2*C_11/C_12]*epsilon = -Pe - Qe = ",
              solcore.asUnit(k.delta_Ehh, 'eV'))
        print("dElh = av * Tre - b[1 + 2*C_11/C_12]*epsilon = -Pe + Qe = ",
              solcore.asUnit(k.delta_Elh, 'eV'))
        print()

    return k
Пример #21
0
def qe_depletion(junction, options):
    """ Calculates the QE curve of a junction object using the depletion approximation as described in J. Nelson, “The Physics of Solar Cells”, Imperial College Press (2003). The junction is then updated with an "iqe" and several "eqe" functions that calculates the QE curve at any wavelength.

    :param junction: A junction object.
    :param options: Solver options.
    :return: None.
    """

    science_reference(
        'Depletion approximation',
        'J. Nelson, “The Physics of Solar Cells”, Imperial College Press (2003).'
    )

    T = options.T

    # First we have to figure out if we are talking about a PN, NP, PIN or NIP junction
    sn = 0 if not hasattr(junction, "sn") else junction.sn
    sp = 0 if not hasattr(junction, "sp") else junction.sp

    # We search for the emitter and check if it is n-type or p-type
    idx = 0
    pn_or_np = 'pn'
    homojunction = True

    for layer in junction:
        if layer.role.lower() != 'emitter':
            idx += 1
        else:
            Na = 0
            Nd = 0
            if hasattr(layer.material, 'Na'): Na = layer.material.Na
            if hasattr(layer.material, 'Nd'): Nd = layer.material.Nd
            if Na < Nd:
                pn_or_np = "np"
                nRegion = junction[idx]
            else:
                pRegion = junction[idx]

            id_top = idx

            break

    # Now we check for an intrinsic region and, if there is, for the base.
    if junction[idx + 1].role.lower() == 'intrinsic':
        iRegion = junction[idx + 1]

        if junction[idx + 2].role.lower() == 'base':
            if pn_or_np == "pn":
                nRegion = junction[idx + 2]
                nidx = idx + 2
            else:
                pRegion = junction[idx + 2]
                pidx = idx + 2

            id_bottom = idx + 2
            homojunction = homojunction and nRegion.material.material_string == pRegion.material.material_string
            homojunction = homojunction and nRegion.material.material_string == iRegion.material.material_string

        else:
            raise RuntimeError(
                'ERROR processing junctions: A layer following the "intrinsic" layer must be defined as '
                '"base".')

    # If there is no intrinsic region, we check directly the base
    elif junction[idx + 1].role.lower() == 'base':
        if pn_or_np == "pn":
            nRegion = junction[idx + 1]
            nidx = idx + 1
        else:
            pRegion = junction[idx + 1]
            pidx = idx + 1
        iRegion = None

        id_bottom = idx + 1
        homojunction = homojunction and nRegion.material.material_string == pRegion.material.material_string

    else:
        raise RuntimeError(
            'ERROR processing junctions: A layer following the "emitter" must be defined as "intrinsic"'
            'or "base".')

    # We assert that we are really working with an homojunction
    assert homojunction, 'ERROR: The DA solver only works with homojunctions, for now.'

    # With all regions identified, it's time to start doing calculations
    kbT = kb * T
    Egap = nRegion.material.band_gap
    R_shunt = min(junction.R_shunt, 1e14) if hasattr(junction,
                                                     'R_shunt') else 1e14

    xp = pRegion.width
    xn = nRegion.width
    xi = 0 if iRegion is None else iRegion.width

    # Now we have to get all the material parameters needed for the calculation
    if hasattr(junction, "dielectric_constant"):
        es = junction.dielectric_constant
    else:
        es = nRegion.material.permittivity  # equal for n and p.  I hope.

    # For the diffusion lenght, subscript n and p refer to the carriers, electrons and holes
    if hasattr(junction, "ln"):
        ln = junction.ln
    else:
        ln = pRegion.material.electron_diffusion_length

    if hasattr(junction, "lp"):
        lp = junction.lp
    else:
        lp = nRegion.material.hole_diffusion_length

    # For the diffusion coefficient, n and p refer to the regions, n side and p side. Yeah, it's confusing...
    if hasattr(junction, "mup"):
        dp = junction.mup * kbT / q
    else:
        dp = pRegion.material.electron_mobility * kbT / q

    if hasattr(junction, "mun"):
        dn = junction.mun * kbT / q
    else:
        dn = nRegion.material.hole_mobility * kbT / q

    # Effective masses and effective density of states
    # mEff_h = nRegion.material.eff_mass_hh_z * electron_mass
    # mEff_e = pRegion.material.eff_mass_electron * electron_mass
    #
    # Nv = 2 * (mEff_h * kb * T / (2 * pi * hbar ** 2)) ** 1.5  # Jenny p58
    # Nc = 2 * (mEff_e * kb * T / (2 * pi * hbar ** 2)) ** 1.5
    # niSquared = Nc * Nv * np.exp(-Egap / (kb * T))
    # ni = np.sqrt(niSquared)

    ni = nRegion.material.ni
    niSquared = ni**2

    Na = pRegion.material.Na
    Nd = nRegion.material.Nd
    Vbi = (kbT / q) * np.log(Nd * Na / niSquared) if not hasattr(
        junction, "Vbi") else junction.Vbi  # Jenny p146

    # It's time to calculate the depletion widths
    if not hasattr(junction, "wp") or not hasattr(junction, "wn"):

        if hasattr(
                junction, "depletion_approximation"
        ) and junction.depletion_approximation == "one-sided abrupt":
            print(
                "using one-sided abrupt junction approximation for depletion width"
            )
            science_reference(
                "Sze abrupt junction approximation",
                "Sze: The Physics of Semiconductor Devices, 2nd edition, John Wiley & Sons, Inc (2007)"
            )
            wp = np.sqrt(2 * es * Vbi / (q * Na))
            wn = np.sqrt(2 * es * Vbi / (q * Nd))

        else:
            wn = (-xi + np.sqrt(xi**2 + 2. * es * Vbi / q *
                                (1 / Na + 1 / Nd))) / (1 + Nd / Na)
            wp = (-xi + np.sqrt(xi**2 + 2. * es * Vbi / q *
                                (1 / Na + 1 / Nd))) / (1 + Na / Nd)

    wn = wn if not hasattr(junction, "wn") else junction.wn
    wp = wp if not hasattr(junction, "wp") else junction.wp

    w = wn + wp + xi

    # Now it is time to calculate currents
    if pn_or_np == "pn":
        l_top, l_bottom = ln, lp
        x_top, x_bottom = xp, xn
        w_top, w_bottom = wp, wn
        s_top, s_bottom = sp, sn
        d_top, d_bottom = dp, dn
        min_top, min_bot = niSquared / Na, niSquared / Nd
    else:
        l_bottom, l_top = ln, lp
        x_bottom, x_top = xp, xn
        w_bottom, w_top = wp, wn
        s_bottom, s_top = sp, sn
        d_bottom, d_top = dp, dn
        min_bot, min_top = niSquared / Na, niSquared / Nd

    widths = []
    for layer in junction:
        widths.append(layer.width)

    cum_widths = np.cumsum([0] + widths)

    g = junction.absorbed
    wl = options.wavelength
    wl_sp, ph = options.light_source.spectrum(output_units='photon_flux_per_m',
                                              x=wl)

    # The contribution from the Emitter (top side).
    xa = cum_widths[id_top]
    xb = cum_widths[id_top + 1] - w_top
    deriv = get_J_sc_diffusion_vs_WL(xa,
                                     xb,
                                     g,
                                     d_top,
                                     l_top,
                                     min_top,
                                     s_top,
                                     wl,
                                     ph,
                                     side='top')
    j_sc_top = d_top * abs(deriv)

    # The contribution from the Base (bottom side).
    xa = cum_widths[id_bottom] + w_bottom
    xb = cum_widths[id_bottom + 1]
    deriv = get_J_sc_diffusion_vs_WL(xa,
                                     xb,
                                     g,
                                     d_bottom,
                                     l_bottom,
                                     min_bot,
                                     s_bottom,
                                     wl,
                                     ph,
                                     side='bottom')
    j_sc_bot = d_bottom * abs(deriv)

    # The contribution from the SCR (includes the intrinsic region, if present).
    xa = cum_widths[id_top + 1] - w_top
    xb = cum_widths[id_bottom] + w_bottom
    j_sc_scr = get_J_sc_SCR_vs_WL(xa, xb, g, wl, ph)

    # The total light absorbed, but not necessarily collected, is:
    xa = cum_widths[id_top]
    xb = cum_widths[id_bottom + 1]
    zz = np.linspace(xa, xb, 10001)
    gg = g(zz) * ph
    current_absorbed = np.trapz(gg, zz, axis=0)

    # Now, we put everything together
    j_sc = j_sc_top + j_sc_bot + j_sc_scr

    eqe = j_sc / ph
    eqe_emitter = j_sc_top / ph
    eqe_base = j_sc_bot / ph
    eqe_scr = j_sc_scr / ph

    junction.iqe = interp1d(wl, j_sc / current_absorbed)
    junction.eqe = interp1d(wl,
                            eqe,
                            kind='linear',
                            bounds_error=False,
                            assume_sorted=True,
                            fill_value=(eqe[0], eqe[-1]))
    junction.eqe_emitter = interp1d(wl,
                                    eqe_emitter,
                                    kind='linear',
                                    bounds_error=False,
                                    assume_sorted=True,
                                    fill_value=(eqe_emitter[0],
                                                eqe_emitter[-1]))
    junction.eqe_base = interp1d(wl,
                                 eqe_base,
                                 kind='linear',
                                 bounds_error=False,
                                 assume_sorted=True,
                                 fill_value=(eqe_base[0], eqe_base[-1]))
    junction.eqe_scr = interp1d(wl,
                                eqe_scr,
                                kind='linear',
                                bounds_error=False,
                                assume_sorted=True,
                                fill_value=(eqe_scr[0], eqe_scr[-1]))

    junction.qe = State({
        'WL': wl,
        'IQE': junction.iqe(wl),
        'EQE': junction.eqe(wl),
        'EQE_emitter': junction.eqe_emitter(wl),
        'EQE_base': junction.eqe_base(wl),
        'EQE_scr': junction.eqe_scr(wl)
    })
Пример #22
0
def test_state():
    array1 = np.array([1, 2, 3])
    array2 = np.array([2, 3, 4])

    state1 = State()
    state1.a = 1
    state1.b = True
    state1.c = array1
    assert state1.__len__() == 3
    assert state1.a == 1
    assert state1.b == True
    assert all([input == output for input, output in zip(array1, state1.c)])

    state2 = State()
    state2.b = False
    state2.d = array2
    state1.update(state2)
    assert state1.__len__() == 4
    assert state1.a == 1
    assert state2.b == False
    assert all([input == output for input, output in zip(array1, state1.c)])
    assert all([input == output for input, output in zip(array2, state1.d)])
def iv_depletion(junction, options):
    """ Calculates the IV curve of a junction object using the depletion approximation as described in J. Nelson, “The Physics of Solar Cells”, Imperial College Press (2003). The junction is then updated with an "iv" function that calculates the IV curve at any voltage.

    :param junction: A junction object.
    :param options: Solver options.
    :return: None.
    """

    science_reference(
        'Depletion approximation',
        'J. Nelson, “The Physics of Solar Cells”, Imperial College Press (2003).'
    )

    junction.voltage = options.internal_voltages
    T = options.T
    kbT = kb * T

    id_top, id_bottom, pRegion, nRegion, iRegion, pn_or_np = identify_layers(
        junction)
    xn, xp, xi, sn, sp, ln, lp, dn, dp, Nd, Na, ni, es = identify_parameters(
        junction, T, pRegion, nRegion, iRegion)

    niSquared = ni**2

    Vbi = (kbT / q) * np.log(Nd * Na / niSquared) if not hasattr(
        junction, "Vbi") else junction.Vbi  # Jenny p146

    #Na, Nd, ni, niSquared, xi, ln, lp, xn, xp, sn, sp, dn, dp, es, id_top, id_bottom, Vbi, pn_or_np = process_junction(junction, options)

    R_shunt = min(junction.R_shunt, 1e14) if hasattr(junction,
                                                     'R_shunt') else 1e14

    # And now we account for the possible applied voltage, which can be, at most, equal to Vbi
    V = np.where(junction.voltage < Vbi - 0.001, junction.voltage, Vbi - 0.001)

    wn, wp = get_depletion_widths(junction, es, Vbi, V, Na, Nd, xi)

    w = wn + wp + xi

    # Now it is time to calculate currents
    if pn_or_np == "pn":
        l_top, l_bottom = ln, lp
        x_top, x_bottom = xp, xn
        w_top, w_bottom = wp, wn
        s_top, s_bottom = sp, sn
        d_top, d_bottom = dp, dn
        min_top, min_bot = niSquared / Na, niSquared / Nd
    else:
        l_bottom, l_top = ln, lp
        x_bottom, x_top = xp, xn
        w_bottom, w_top = wp, wn
        s_bottom, s_top = sp, sn
        d_bottom, d_top = dp, dn
        min_bot, min_top = niSquared / Na, niSquared / Nd

    JtopDark = get_j_dark(x_top, w_top, l_top, s_top, d_top, V, min_top, T)
    JbotDark = get_j_dark(x_bottom, w_bottom, l_bottom, s_bottom, d_bottom, V,
                          min_bot, T)

    # hereby we define the subscripts to refer to the layer in which the current is generated:
    if pn_or_np == "pn":
        JnDark, JpDark = JbotDark, JtopDark
    else:
        JpDark, JnDark = JbotDark, JtopDark

    # These might not be the right lifetimes. Actually, they are not as they include all recombination processes, not
    # just SRH recombination, which is what the equation in Jenny, p159 refers to. Let´ leave them, for now.
    lifetime_n = ln**2 / dn
    lifetime_p = lp**2 / dp  # Jenny p163

    # Here we use the full version of the SRH recombination term as calculated by Sah et al. Works for positive bias
    # and moderately negative ones.

    science_reference(
        'SRH current term.',
        'C. T. Sah, R. N. Noyce, and W. Shockley, “Carrier Generation and Recombination in P-N Junctions and P-N Junction Characteristics,” presented at the Proceedings of the IRE, 1957, vol. 45, no. 9, pp. 1228–1243.'
    )
    Jrec = get_Jsrh(ni, V, Vbi, lifetime_p, lifetime_n, w, kbT)

    J_sc_top = 0
    J_sc_bot = 0
    J_sc_scr = 0

    if options.light_iv:

        widths = []
        for layer in junction:
            widths.append(layer.width)

        cum_widths = np.cumsum([0] + widths)

        g = junction.absorbed
        wl = options.wavelength
        wl_sp, ph = options.light_source.spectrum(
            output_units='photon_flux_per_m', x=wl)
        id_v0 = np.argmin(abs(V))

        # The contribution from the Emitter (top side).
        xa = cum_widths[id_top]
        xb = cum_widths[id_top + 1] - w_top[id_v0]

        if options.da_mode == 'bvp':
            deriv = get_J_sc_diffusion(xa,
                                       xb,
                                       g,
                                       d_top,
                                       l_top,
                                       min_top,
                                       s_top,
                                       wl,
                                       ph,
                                       side='top')
        else:
            xbb = xb - (xb - xa) / 1001.
            deriv = get_J_sc_diffusion_green(xa,
                                             xbb,
                                             g,
                                             d_top,
                                             l_top,
                                             min_top,
                                             s_top,
                                             ph,
                                             side='top')
            deriv = np.trapz(deriv, wl)
        J_sc_top = q * d_top * abs(deriv)

        # The contribution from the Base (bottom side).
        xa = cum_widths[id_bottom] + w_bottom[id_v0]
        xb = cum_widths[id_bottom + 1]
        if options.da_mode == 'bvp':
            deriv = get_J_sc_diffusion(xa,
                                       xb,
                                       g,
                                       d_bottom,
                                       l_bottom,
                                       min_bot,
                                       s_bottom,
                                       wl,
                                       ph,
                                       side='bottom')
        else:
            xbb = xb - (xb - xa) / 1001.
            deriv = get_J_sc_diffusion_green(xa,
                                             xbb,
                                             g,
                                             d_bottom,
                                             l_bottom,
                                             min_bot,
                                             s_bottom,
                                             ph,
                                             side='bottom')
            deriv = np.trapz(deriv, wl)
        J_sc_bot = q * d_bottom * abs(deriv)

        # The contribution from the SCR (includes the intrinsic region, if present).
        xa = cum_widths[id_top + 1] - w_top[id_v0]
        xb = cum_widths[id_bottom] + w_bottom[id_v0]
        J_sc_scr = q * get_J_sc_SCR(xa, xb, g, wl, ph)

    # And, finally, we output the currents
    junction.current = Jrec + JnDark + JpDark + V / R_shunt - J_sc_top - J_sc_bot - J_sc_scr
    junction.iv = interp1d(junction.voltage,
                           junction.current,
                           kind='linear',
                           bounds_error=False,
                           assume_sorted=True,
                           fill_value=(junction.current[0],
                                       junction.current[-1]))
    junction.region_currents = State({
        "Jn_dif": JnDark,
        "Jp_dif": JpDark,
        "Jscr_srh": Jrec,
        "J_sc_top": J_sc_top,
        "J_sc_bot": J_sc_bot,
        "J_sc_scr": J_sc_scr
    })
Пример #24
0
from solcore.structure import Layer, Junction, TunnelJunction
from solcore.absorption_calculator import calculate_rat_rcwa, calculate_absorption_profile_rcwa

import numpy as np
import types
from solcore.state import State


rcwa_options = State()
rcwa_options.size = [500, 500]
rcwa_options.orders = 4
rcwa_options.theta = 0
rcwa_options.phi = 0
rcwa_options.pol = 'u'
rcwa_options.parallel = False
rcwa_options.n_jobs = -1

def solve_rcwa(solar_cell, options):
    """Calculates the reflection, transmission and absorption of a solar cell object using the rigorous coupled-wave analysis solver.

    :param solar_cell: A solar_cell object
    :param options: Options for the solver
    :return: None
    """
    wl = options.wavelength
    solar_cell.wavelength = options.wavelength
    parallel = options.parallel
    rcwa_options = options.get("rcwa_options", None)

    # We include the shadowing losses
    initial = (1 - solar_cell.shading) if hasattr(solar_cell, 'shading') else 1
Пример #25
0
import solcore.analytic_solar_cells as ASC
from solcore.light_source import LightSource
from solcore.state import State
from solcore.optics import solve_beer_lambert, solve_tmm, solve_rcwa, rcwa_options, solve_external_optics
from solcore.absorption_calculator import RCWASolverError
from solcore.structure import Layer, Junction, TunnelJunction

try:
    import solcore.poisson_drift_diffusion as PDD

    a = PDD.pdd_options
except AttributeError:
    PDD.pdd_options = {}

default_options = State()
pdd_options = PDD.pdd_options
asc_options = ASC.db_options


def merge_dicts(*dict_args):
    """
    Given any number of dicts, shallow copy and merge into a new dict,
    precedence goes to key value pairs in latter dicts.
    """
    result = State()
    for dictionary in dict_args:
        result.update(dictionary)
    return result

def qe_depletion(junction, options):
    """ Calculates the QE curve of a junction object using the depletion approximation as described in J. Nelson, “The Physics of Solar Cells”, Imperial College Press (2003). The junction is then updated with an "iqe" and several "eqe" functions that calculates the QE curve at any wavelength.

    :param junction: A junction object.
    :param options: Solver options.
    :return: None.
    """

    science_reference(
        'Depletion approximation',
        'J. Nelson, “The Physics of Solar Cells”, Imperial College Press (2003).'
    )

    # First we have to figure out if we are talking about a PN, NP, PIN or NIP junction
    T = options.T
    kbT = kb * T

    id_top, id_bottom, pRegion, nRegion, iRegion, pn_or_np = identify_layers(
        junction)
    xn, xp, xi, sn, sp, ln, lp, dn, dp, Nd, Na, ni, es = identify_parameters(
        junction, T, pRegion, nRegion, iRegion)

    niSquared = ni**2

    Vbi = (kbT / q) * np.log(Nd * Na / niSquared) if not hasattr(
        junction, "Vbi") else junction.Vbi  # Jenny p146

    wn, wp = get_depletion_widths(junction, es, Vbi, 0, Na, Nd, xi)

    # Now it is time to calculate currents
    if pn_or_np == "pn":
        l_top, l_bottom = ln, lp
        w_top, w_bottom = wp, wn
        s_top, s_bottom = sp, sn
        d_top, d_bottom = dp, dn
        min_top, min_bot = niSquared / Na, niSquared / Nd
    else:
        l_bottom, l_top = ln, lp
        w_bottom, w_top = wp, wn
        s_bottom, s_top = sp, sn
        d_bottom, d_top = dp, dn
        min_bot, min_top = niSquared / Na, niSquared / Nd

    widths = []
    for layer in junction:
        widths.append(layer.width)

    cum_widths = np.cumsum([0] + widths)

    g = junction.absorbed
    wl = options.wavelength
    wl_sp, ph = LightSource(source_type='black body', x=wl,
                            T=6000).spectrum(output_units='photon_flux_per_m',
                                             x=wl)

    # The contribution from the Emitter (top side).
    xa = cum_widths[id_top]
    xb = cum_widths[id_top + 1] - w_top

    if options.da_mode == 'bvp':
        deriv = get_J_sc_diffusion_vs_WL(xa,
                                         xb,
                                         g,
                                         d_top,
                                         l_top,
                                         min_top,
                                         s_top,
                                         wl,
                                         ph,
                                         side='top')
    else:
        xbb = xb - (xb - xa) / 1001.
        deriv = get_J_sc_diffusion_green(xa,
                                         xbb,
                                         g,
                                         d_top,
                                         l_top,
                                         min_top,
                                         s_top,
                                         ph,
                                         side='top')
    j_sc_top = d_top * abs(deriv)

    # The contribution from the Base (bottom side).
    xa = cum_widths[id_bottom] + w_bottom
    xb = cum_widths[id_bottom + 1]

    if options.da_mode == 'bvp':
        deriv = get_J_sc_diffusion_vs_WL(xa,
                                         xb,
                                         g,
                                         d_bottom,
                                         l_bottom,
                                         min_bot,
                                         s_bottom,
                                         wl,
                                         ph,
                                         side='bottom')
    else:
        xbb = xb - (xb - xa) / 1001.
        deriv = get_J_sc_diffusion_green(xa,
                                         xbb,
                                         g,
                                         d_bottom,
                                         l_bottom,
                                         min_bot,
                                         s_bottom,
                                         ph,
                                         side='bottom')
    j_sc_bot = d_bottom * abs(deriv)

    # The contribution from the SCR (includes the intrinsic region, if present).
    xa = cum_widths[id_top + 1] - w_top
    xb = cum_widths[id_bottom] + w_bottom
    j_sc_scr = get_J_sc_SCR_vs_WL(xa, xb, g, wl, ph)

    # The total light absorbed, but not necessarily collected, is:
    xa = cum_widths[id_top]
    xb = cum_widths[id_bottom + 1]
    zz = np.linspace(xa, xb, 10001)
    gg = g(zz) * ph
    current_absorbed = np.trapz(gg, zz, axis=0)

    # Now, we put everything together
    j_sc = j_sc_top + j_sc_bot + j_sc_scr

    eqe = j_sc / ph
    eqe_emitter = j_sc_top / ph
    eqe_base = j_sc_bot / ph
    eqe_scr = j_sc_scr / ph

    iqe = j_sc / current_absorbed
    iqe[np.isnan(
        iqe
    )] = 0  # if zero current_absorbed, get NaN in previous line; want 0 IQE

    junction.iqe = interp1d(wl, iqe)

    junction.eqe = interp1d(wl,
                            eqe,
                            kind='linear',
                            bounds_error=False,
                            assume_sorted=True,
                            fill_value=(eqe[0], eqe[-1]))
    junction.eqe_emitter = interp1d(wl,
                                    eqe_emitter,
                                    kind='linear',
                                    bounds_error=False,
                                    assume_sorted=True,
                                    fill_value=(eqe_emitter[0],
                                                eqe_emitter[-1]))
    junction.eqe_base = interp1d(wl,
                                 eqe_base,
                                 kind='linear',
                                 bounds_error=False,
                                 assume_sorted=True,
                                 fill_value=(eqe_base[0], eqe_base[-1]))
    junction.eqe_scr = interp1d(wl,
                                eqe_scr,
                                kind='linear',
                                bounds_error=False,
                                assume_sorted=True,
                                fill_value=(eqe_scr[0], eqe_scr[-1]))

    junction.qe = State({
        'WL': wl,
        'IQE': junction.iqe(wl),
        'EQE': junction.eqe(wl),
        'EQE_emitter': junction.eqe_emitter(wl),
        'EQE_base': junction.eqe_base(wl),
        'EQE_scr': junction.eqe_scr(wl)
    })
import math
from scipy.interpolate import interp1d
import numpy as np

from solcore.constants import kb, q, c, h
from solcore import si
from solcore.state import State

db_options = State()
db_options.db_mode = 'boltzmann'


def iv_detailed_balance(junction, options):
    """ Calculates the IV curve of a junction of kind "DB". This is a detailed balanced calculation of the IV curve of a PN junction characterized by a certain eqe, temperature and chemical potential (voltage). Normally, the eqe will be calculated based on a a given absorption edge and absorptance level resulting in a "top hat" type of eqe, but it can also be provided externally. By default, the solver uses the Boltzmann aproximation, although the full Planck equation might be used, also.

    :param junction: A Junction object of kind "DB"
    :param options: Other arguments for the calculation of the IV.
    :return:
    """
    T = options.T
    light = options.light_iv
    mode = options.db_mode
    Ta = options.T_ambient
    wl = options.wavelength

    try:
        Eg = junction.Eg
        n = junction.n

        if not hasattr(junction, 'eqe'):
            qe_detailed_balance(junction, wl)
Пример #28
0
        Layer(si("200nm"), material=bot_cell_n_material, role='emitter'),
        Layer(si("29800nm"), material=bot_cell_p_material, role='base')
    ],
             sn=1,
             sp=1,
             kind='DA')
]
# And, finally, we put everything together, adding also the surface recombination velocities sn and sp.
# setting kind = 'DA' in the Junction definition tells the electrical solver later to use the depletion approximation
optical_struct = SolarCell(ARC + top_junction + middle_junction + DBRa + DBRb +
                           DBRc + bottom_junction,
                           shading=0.05)

wl = np.linspace(250, 1700, 400) * 1e-9

options = State()
options.wavelength = wl
options.optics_method = 'TMM'
options.no_back_reflection = False
options.pol = 'p'
options.BL_correction = True
options.coherency_list = 111 * ['c']
options.theta = 30
solar_cell_solver(optical_struct, 'qe', options)

plt.figure()
plt.plot(
    wl * 1e9,
    optical_struct[0].layer_absorption + optical_struct[1].layer_absorption)
plt.plot(wl * 1e9, optical_struct[2].layer_absorption)
plt.plot(wl * 1e9, optical_struct[3].layer_absorption)
Пример #29
0
optical_struct = SolarCell([
    Layer(material=GaInP, width=si("5000nm")),
    Junction(
        [
            Layer(material=GaAs, width=si("200nm")),
            Layer(material=GaAs, width=si("5um")),
        ],
        kind="DA",
    ),
    Layer(material=Ge, width=si("50um")),
])

wl = np.linspace(250, 1700, 300) * 1e-9

options = State()
options.wavelength = wl
options.optics_method = "TMM"
options.no_back_reflection = False
options.BL_correction = True
options.recalculate_absorption = True
options.positions = [1e-8, 1e-9, 1e-8, 1e-7]
options.theta = 0

c_list = [
    ["c", "c", "c", "c"],
    ["c", "c", "c", "i"],
    ["c", "i", "i", "c"],
    ["i", "i", "i", "i"],
]
Пример #30
0
def iv_multijunction(solar_cell, options):
    """ Calculates the overall IV characteristics of any number of junctions numerically at the requested voltage
    points. If photocurrent is not provided, the resultant IV characteristics are purely recombination currents,
    otherwise light IVs are returned.

    In the end, the soalr_cell object is updated with an "iv" attribute containing a dictionary with:
        "IV": (V, I) Calculated IV characteristics
        "junction IV": [(V junc 1, I junc 1), (V junc 2, I junc 2), ...]
        "Rseries IV": (V, I) Calculated IV characteristics of the series resistance
        "coupled IV": For each junction but for the first one, a list with the coupled current coming form the upper junctions.
        "Isc", Voc", "P", "FF" and "Eta": In case of mpp = True and light IV.

    The sign convention is:
        - Photocurrents: Positive.
        - Dark Currents: Negative

    :param solar_cell: A solar cell object with one or more junctions. The IV of the individual junctions must have been calculated already.
    :param kwargs: A dictionary containing, at least, the following elements:
        - mpp: (Boolean) If Isc, Voc, FF, Vmpp, Impp and Pmpp must be calculated.
        - voltages: Array of voltages in which to calculate the data
        - light_iv: (Boolean) if light IV is being calculated
    :return: None

    """
    output_V = options.voltages
    mpp_parameters = options.mpp
    radiative_coupling = options.radiative_coupling
    light_iv = options.light_iv

    Rs = max(solar_cell.R_series, 1e-14)

    # The current and junction voltage arrays
    num_jun = solar_cell.junctions
    tunnel_jun = solar_cell.tunnel_indices
    V_junction_array = np.zeros((len(output_V), num_jun))

    # The following assumes that all junctions have the currents defined at the same voltages
    minimum_J = solar_cell(0).current
    for j in range(1, num_jun):
        minimum_J = np.where(
            np.absolute(minimum_J) < np.absolute(solar_cell(j).current),
            minimum_J,
            solar_cell(j).current)

    minimum_J = np.sort(minimum_J)

    # For each junction, we find the voltage corresponding to each of these currents
    temp_V_junction_array = np.zeros((len(minimum_J), num_jun))
    for j in range(num_jun):
        temp_V_junction_array[:, j] = np.interp(minimum_J,
                                                solar_cell(j).current,
                                                solar_cell(j).voltage)

    # We calculate the total voltage related to the...
    # ... series resistance
    temp_V_total = Rs * minimum_J

    # ... tunnel junctions
    for j in tunnel_jun:
        temp_V_total -= solar_cell[j].vi(-minimum_J)

    # ... and the normal junctions
    for j in range(num_jun):
        temp_V_total += temp_V_junction_array[:, j]

    # Finally, we calculate the current at the requested voltages and the voltages for each junction
    output_J = np.interp(output_V, temp_V_total, minimum_J)
    for j in range(num_jun):
        V_junction_array[:, j] = np.interp(output_V, temp_V_total,
                                           temp_V_junction_array[:, j])

    coupled_J = []
    if radiative_coupling:
        # First of all, we check if all the junctions are DB or 2D
        ok = True
        for i in range(num_jun):
            ok = ok and solar_cell(i).kind in ['DB', '2D']

        if ok:
            output_J, V_junction_array, coupled_J = solve_radiative_coupling(
                solar_cell, options, V_junction_array)
        else:
            print(
                'WARNING: Only "DB" and "2D" junctions (using the DB solver) can use radiative coupling.\nSkipping calculation...'
            )

    # Finally, we calculate the solar cell parameters
    Isc = None
    Voc = None
    Pmpp = None
    Vmpp = None
    Impp = None
    FF = None
    Eta = None

    # If we are calculating the light IV, we also calculate the main parameters: Jsc, Voc, FF, MPP...
    if mpp_parameters and light_iv:
        VV = output_V
        II = -output_J
        PP = VV * II
        power = options.light_source.power_density
        try:
            Isc = np.interp([0], VV, II)[0]
            Voc = np.interp([0], -II, VV)[0]
            Pmpp = np.max(PP)
            idx = np.argmax(PP)
            Vmpp = VV[idx]
            Impp = II[idx]
            FF = Pmpp / (Isc * Voc)
            Eta = Pmpp / power

        except Exception as err:
            print('Error calculating the MPP parameters: {}'.format(err))

    solar_cell.iv = State({
        "IV": (output_V, -output_J),
        "junction IV":
        [(V_junction_array[:, i], -output_J) for i in range(num_jun)],
        "coupled IV":
        coupled_J,
        "Rseries IV": (output_J * Rs, -output_J),
        "Isc":
        Isc,
        "Voc":
        Voc,
        "FF":
        FF,
        "Pmpp":
        Pmpp,
        "Vmpp":
        Vmpp,
        "Impp":
        Impp,
        "Eta":
        Eta
    })