pySIMsalabim.experiments package
Submodules
pySIMsalabim.experiments.CV module
Perform capacitance simulations
- pySIMsalabim.experiments.CV.MottSchottky_plot(session_path, output_file, xscale='linear', yscale='linear', plot_type=<function errorbar>)[source]
Plot the 1/capacitance^2 against voltage
- Parameters:
session_path (string) – working directory for zimt
output_file (string) – Filename where the capacitance data is stored
xscale (string) – Scale of the x-axis. E.g linear or log
yscale (string) – Scale of the y-axis. E.g linear or log
plot_type (matplotlib.pyplot) – Type of plot to display
- pySIMsalabim.experiments.CV.calc_capacitance(data, del_V, freq)[source]
Calculate the capacitance over the time range
- Parameters:
data (dataFrame) – Pandas dataFrame containing the time, voltage, current density and numerical error in the current density of the tj_file
del_V (float) – Voltage step that is applied directly after t=0
freq (float) – Frequency at which the capacitance-voltage measurement is performed
- Returns:
np.array – Array of the capacitance
np.array – Array of the capacitance error
- pySIMsalabim.experiments.CV.calc_capacitance_forOneVoltage(I, errI, time, VStep, freq)[source]
Fourier Decomposition formula which computes the capacitance at frequency freq (Hz) and its complex error Based on S.E. Laux, IEEE Trans. Electron Dev. 32 (10), 2028 (1985), eq. 5b We integrate from 0 to time[imax] at frequency 1/time[imax]
- Parameters:
I (np.array) – Array of currents
errI (np.array) – Numerical error in calculated currents (output of ZimT)
time (np.array) – Array with all time positions, from 0 to tmax
VStep (float) – Voltage step
imax (integer) – Index of the last timestep/first frequency for which the integrals are calculated
- Returns:
float – Capacitance at frequency f: C(f)
float – Numerical error in calculated capacitance
- pySIMsalabim.experiments.CV.calc_impedance_CV(data, del_V, isToPlot, session_path, zimt_device_parameters, Rseries=0, Rshunt=-1000.0)[source]
Calculate the impedance over the frequency range
- Parameters:
data (dataFrame) – Pandas dataFrame containing the time, voltage, current density and numerical error in the current density of the tj_file
del_V (float) – Voltage step
isToPlot (list) – List of array indices that will be used in the plotting
- Returns:
np.array – Array of frequencies
np.array – Array of the real component of impedance
np.array – Array of the imaginary component of impedance
np.array – Array of complex error
np.array – Array of capacitance
np.array – Array of conductance
np.array – Array of error in capacitance
np.array – Array of error in conductance
- pySIMsalabim.experiments.CV.calc_impedance_limit_time(I, errI, time, VStep, imax)[source]
Fourier Decomposition formula which computes the impedance at frequency freq (Hz) and its complex error Based on S.E. Laux, IEEE Trans. Electron Dev. 32 (10), 2028 (1985), eq. 5a, 5b We integrate from 0 to time[imax] at frequency 1/time[imax]
- Parameters:
I (np.array) – Array of currents
errI (np.array) – Numerical error in calculated currents (output of ZimT)
time (np.array) – Array with all time positions, from 0 to tmax
VStep (float) – Voltage step
imax (integer) – Index of the last timestep/first frequency for which the integrals are calculated
- Returns:
float – Frequency belonging to Z(f)
complex number – Impedance at frequency f: Z(f)
complex number – Numerical error in calculated impedance Z(f)
- pySIMsalabim.experiments.CV.cap_plot(session_path, output_file, xscale='linear', yscale='linear', plot_type=<function errorbar>)[source]
Plot the capacitance against voltage
- Parameters:
session_path (string) – working directory for zimt
output_file (string) – Filename where the capacitance data is stored
xscale (string) – Scale of the x-axis. E.g linear or log
yscale (string) – Scale of the y-axis. E.g linear or log
plot_type (matplotlib.pyplot) – Type of plot to display
- pySIMsalabim.experiments.CV.create_tVG_CV(V_0, V_max, del_V, V_step, G_frac, tVG_name, session_path, freq, ini_timeFactor, timeFactor)[source]
Create a tVG file for capacitance experiment.
- Parameters:
V_0 (float) – Initial voltage
V_max (float) – Maximum voltage
del_V (float) – Voltage step that is applied directly after t=0
V_step (float) – Voltage difference, determines at which voltages the capacitance is determined
G_frac (float) – Fractional light intensity
tVG_name (string) – Name of the tVG file
session_path (string) – Path of the simulation folder for this session
freq (float) – Frequency at which the capacitance-voltage measurement is performed
ini_timeFactor (float) – Constant defining the size of the initial timestep
timeFactor (float) – Exponential increase of the timestep, to reduce the amount of timepoints necessary. Use values close to 1.
- Returns:
A message to indicate the result of the process
- Return type:
string
- pySIMsalabim.experiments.CV.create_tVG_CV_Rseries(Vint, del_V, G_frac, tVG_name, session_path, freq, ini_timeFactor, timeFactor)[source]
Create a tVG file for capacitance experiment.
- Parameters:
V_int (float) – Adjusted internal voltage to correct for series resistance later
del_V (float) – Voltage step that is applied directly after t=0
G_frac (float) – Fractional light intensity
tVG_name (string) – Name of the tVG file
session_path (string) – Path of the simulation folder for this session
freq (float) – Frequency at which the capacitance-voltage measurement is performed
ini_timeFactor (float) – Constant defining the size of the initial timestep
timeFactor (float) – Exponential increase of the timestep, to reduce the amount of timepoints necessary. Use values close to 1.
- Returns:
A message to indicate the result of the process
- Return type:
string
- pySIMsalabim.experiments.CV.create_tVG_SS_CV(V_min, V_max, Vstep, G_frac, tVG_name, session_path)[source]
Creates the tVG file for the steady state simulation with only t 0 and V_0
- Parameters:
V_min (float) – Initial voltage
V_max (float) – Maximum voltage
del_V (float) – Voltage step that is applied after t=0
G_frac (float) – Fractional light intensity
tVG_name (string) – Name of the tVG file
session_path (string) – Path of the simulation folder for this session
- Returns:
A message to indicate the result of the process
- Return type:
string
- pySIMsalabim.experiments.CV.get_capacitance(data, freq, V_0, V_max, del_V, V_step, session_path, zimt_device_parameters, output_file)[source]
Calculate the capacitance from the simulation result
- Parameters:
data (DataFrame) – DataFrame with the simulation results (tj.dat) file
freq (float) – Frequency at which the capacitance-voltage measurement is performed
V_0 (float) – Initial voltage
V_max (float) – Maximum voltage
del_V (float) – Voltage step that is applied directly after t=0
V_step (float) – Voltage difference, determines at which voltages the capacitance is determined
session_path (string) – working directory for zimt
output_file (string) – Filename where the capacitance data is stored
- Returns:
returns -1 (failed) or 1 (success), including a message
- Return type:
integer,string
- pySIMsalabim.experiments.CV.plot_capacitance(session_path, output_file='CapVol.dat')[source]
Make a plot of the capacitance against voltage and the Mott-Schottky plot
- Parameters:
session_path (string) – working directory for zimt
- pySIMsalabim.experiments.CV.run_CV_simu(zimt_device_parameters, session_path, freq, V_min, V_max, V_step, G_frac=1, del_V=0.01, run_mode=False, tVG_name='tVG.txt', output_file='CapVol.dat', tj_name='tj.dat', varFile='none', ini_timeFactor=0.001, timeFactor=1.02, **kwargs)[source]
Create a tVG file and run ZimT with capacitance device parameters
- Parameters:
zimt_device_parameters (string) – Name of the zimt device parameters file
session_path (string, optional) – working directory for zimt
freq (float) – Frequency at which the capacitance-voltage measurement is performed
V_min (float) – Initial voltage
V_max (float) – Maximum voltage
V_step (float) – Voltage difference, determines at which voltages the capacitance is determined
G_frac (float, optional) – Fractional light intensity, by default 1
del_V (float, optional) – Voltage step that is applied directly after t=0, by default 0.01
run_mode (bool, optional) – Indicate whether the script is in ‘web’ mode (True) or standalone mode (False). Used to control the console output, by default False
tVG_name (string, optional) – Name of the tVG file, by default tVG.txt
output_file (string, optional) – Name of the file where the capacitance data is stored, by default CapVol.dat
tj_name (string, optional) – Name of the tj file where the capacitance data is stored, by default tj.dat
varFile (string, optional) – Name of the var file, by default ‘none’
ini_timeFactor (float, optional) – Constant defining the size of the initial timestep, by default 1e-3
timeFactor (float, optional) – Exponential increase of the timestep, to reduce the amount of timepoints necessary. Use values close to 1., by default 1.02
- Returns:
CompletedProcess – Output object of with returncode and console output of the simulation
string – Return message to display on the UI, for both success and failed
- pySIMsalabim.experiments.CV.store_capacitance_data(session_path, V, cap, errC, output_file='CapVol.dat')[source]
Save the capacitance & its error in one file called CapVol.dat
- Parameters:
session_path (string) – working directory for zimt
V (np.array) – Array of frequencies
cap (np.array) – Array of the capacitance
errC (np.array) – Array of the capacitance error
output_file (string) – Filename where the capacitance data is stored
pySIMsalabim.experiments.EQE module
Perform EQE simulations
- pySIMsalabim.experiments.EQE.EQE_create_spectrum_files(lambda_array, p, session_path, spectrum_path, tmp_spectrum_path)[source]
Create all the spectrum files with monochromatic peaks at the specified wavelength, with a fixed height set by the number of added photons p. The files are stored in a temporary folder: tmp_spectrum.
- Parameters:
lambda_array (array) – Array with the wavelengths at which the monochromatic peaks will be created.
p (float) – Number of photons that are added to the irradiance.
session_path (string) – Path to the session folder where the files must be created.
spectrum_path (string) – Path to the original spectrum file.
tmp_spectrum_path (string) – Path to the temporary folder where the modified spectrum files will be stored.
- pySIMsalabim.experiments.EQE.calc_EQE(Jext, Jext_err, J0_single, J0_err_single, lambda_array, p)[source]
Calculate the EQE values for the monochromatic peaks at each wavelength. Based on the change in the short-circuit current and the number of added photons.
- Parameters:
Jext (array) – Short-circuit current density for each monochromatic peak.
Jext_err (array) – Error in the short-circuit current density for each monochromatic peak.
J0_single (float) – Short-circuit current density for the normal spectrum.
J0_err_single (float) – Error in the short-circuit current density for the normal spectrum.
lambda_array (array) – Array with the wavelengths at which the monochromatic peaks were created.
p (float) – Number of photons that are added to the irradiance.
- Returns:
Arrays with the change in short-circuit current density, error in the change in short-circuit currentdensity, monochromatic intensity, EQE values and error in the EQE values.
- Return type:
array, array, array, array, array
- pySIMsalabim.experiments.EQE.get_CurrDens(JV_file, session_path)[source]
Get the current density and its from the JV_file as stored in the first row. :param JV_file: Name of the file where the JV data is stored. Must be unique for each simulation, by default ‘JV.dat’ :type JV_file: str, optional :param session_path: Path to the session folder where the simulation will run. :type session_path: string
- Returns:
Short-circuit current and its error.
- Return type:
float, float
- pySIMsalabim.experiments.EQE.run_EQE(simss_device_parameters, session_path, spectrum, lambda_min, lambda_max, lambda_step, Vext, output_file='EQE.dat', JV_file_name='JV.dat', varFile='none', remove_dirs=True, parallel=False, max_jobs=3, run_mode=True, **kwargs)[source]
Run the EQE calculation for a given spectrum and external voltage, and save the results in a file.
- Parameters:
simss_device_parameters (string) – Name of the device parameters file.
session_path (string) – Path to the session folder where the simulation will run.
spectrum (string) – Path to the original spectrum file.
lambda_min (float) – Minimum wavelength for the spectrum.
lambda_max (float) – Maximum wavelength for the spectrum.
lambda_step (float) – Step size for the wavelength.
Vext (float) – External voltage at which the simulation will run and the EQE must be calculated.
output_file (string, optional) – Name of the file where the results will be stored. The file will be stored in the session folder, by default ‘EQE.dat’
JV_file_name (string, optional) – Name of the JV file. Must be unique for each simulation, by default ‘JV.dat’
varFile (string, optional) – Name of the var file, by default ‘none’
remove_dirs (bool, optional) – Remove the temporary directories, by default True
parallel (bool, optional) – Run the simulations in parallel, by default False
max_jobs (int, optional) – Maximum number of parallel jobs, by default max(1,os.cpu_count()-1)
run_mode (bool, optional) – indicate whether the script is in ‘web’ mode (True) or standalone mode (False). Used to control the console output, by default True
**kwargs (dict) – Additional arguments to be passed to the function.
- Returns:
0 if the function runs successfully.
- Return type:
int
pySIMsalabim.experiments.JV_steady_state module
Perform steady-state JV simulations
- pySIMsalabim.experiments.JV_steady_state.run_SS_JV(simss_device_parameters, session_path, JV_file_name='JV.dat', varFile='none', G_fracs=[], parallel=False, max_jobs=3, run_mode=True, **kwargs)[source]
- Parameters:
simss_device_parameters (string) – Name of the device parameters file.
session_path (string) – Path to the session folder where the simulation will run.
JV_file_name (string) – Name of the JV file.
parallel (bool, optional) – Run the simulations in parallel, by default False
max_jobs (int, optional) – Maximum number of parallel jobs, by default max(1,os.cpu_count()-1)
cmd_pars (_type_, optional) – _description_, by default None
UUID (str, optional) – _description_, by default ‘’
run_mode (bool, optional) – indicate whether the script is in ‘web’ mode (True) or standalone mode (False). Used to control the console output, by default True
**kwargs (dict) – Additional arguments to be passed to the function.
- Returns:
int – Exitcode of the simulation.
str – Message from the simulation.
pySIMsalabim.experiments.hysteresis module
Perform JV hysteresis simulations
- pySIMsalabim.experiments.hysteresis.Compare_Exp_Sim_JV(session_path, expJV_Vmin_Vmax, expJV_Vmax_Vmin, rms_mode, direction, tj_file_name='tj.dat')[source]
Calculate the root-mean-square (rms) error of the simulated data compared to the experimental data. The used formulas are described in the Manual (see the variable rms_mode in the section ‘Description of simulated device parameters’).
- Parameters:
session_path (string) – Path of the simulation folder for this session
expJV_Vmin_Vmax (string) – Name of the file of the forward JV scan
expJV_Vmax_Vmin (string) – Name of the file of the backward JV scan
rms_mode (string) – Indicates how the normalised rms error should be calculated: either in linear or logarithmic form
direction (integer) – Perform a Vmin-Vmax-Vmin (1) or Vmax-Vmin-Vmax scan (-1)
tj_file_name (dataFrame) – Pandas dataFrame containing the tj output file from ZimT
- Returns:
Calculated rms-error
- Return type:
Float
- pySIMsalabim.experiments.hysteresis.Hysteresis_JV(zimt_device_parameters, session_path, UseExpData, scan_speed, direction, G_frac, tVG_name='tVG.txt', tj_name='tj.dat', varFile='none', run_mode=False, Vmin=0.0, Vmax=0.0, steps=0, expJV_Vmin_Vmax='', expJV_Vmax_Vmin='', rms_mode='lin', **kwargs)[source]
Create a tVG file and perform a JV hysteresis experiment.
- Parameters:
zimt_device_parameters (string) – name of the zimt device parameters file
session_path (string) – working directory for zimt
UseExpData (integer) – If 1, use experimental JV curves. If 0, Use Vmin, Vmax as boundaries
scan_speed (float) – Voltage scan speed [V/s]
direction (integer) – Perform a forward-backward (1) or backward-forward scan (-1).
G_frac (float) – Device Parameter | Fractional generation rate
run_mode (bool, optional) – Indicate whether the script is in ‘web’ mode (True) or standalone mode (False). Used to control the console output, by default False
tVG_name (string, optional) – Device Parameter | Name of the tVG file, by default ‘tVG.txt’
tj_name (string, optional) – Name of the tj file, by default ‘tj.dat’
varFile (string, optional) – Name of the var file, by default ‘none’
run_mode – indicate whether the script is in ‘web’ mode (True) or standalone mode (False). Used to control the console output, by default False
Vmin (float, optional) – Left voltage boundary, by default 0.0
Vmax (float, optional) – Right voltage boundary, by default 0.0
steps (int, optional) – Number of time steps, by default 0
expJV_Vmin_Vmax (str, optional) – file name of the first expJV curve, by default ‘’
expJV_Vmax_Vmin (str, optional) – file name of the second expJV curve, by default ‘’
rms_mode (str, optional) – Either ‘lin’ or ‘log’ to specify how the rms error is calculated
- Returns:
CompletedProcess – Output object of with returncode and console output of the simulation
string – Return message to display on the UI, for both success and failed
dict – Dictionary containing the special output values of the simulation. In this case, the rms error (‘rms’) and the hysteresis index (‘hyst_index’)
- pySIMsalabim.experiments.hysteresis.build_tVG_arrays(Vmin, Vmax, scan_speed, direction, steps, G_frac)[source]
Build the Arrays for time, voltage and Generation rate for a hysteresis experiment.
- Parameters:
Vmin (float) – minimum voltage
Vmax (float) – maximum voltage
scan_speed (float) – Voltage scan speed [V/s]
direction (integer) – Perform a Vmin-Vmax-Vmin (1) or Vmax-Vmin-Vmax scan (-1)
steps (integer) – Number of time steps
G_frac (float) – Device Parameter | Fractional Generation rate
- Returns:
np.array – Array of time points
np.array – Array of voltages
np.array – Array of generation rates
- pySIMsalabim.experiments.hysteresis.build_tVG_arrays_log(Vmin, Vmax, scan_speed, direction, steps, G_frac, Vminexpo=0.01)[source]
Build the Arrays for time, voltage and Generation rate for a hysteresis experiment with an exponential voltage sweep.
- Parameters:
Vmin (float) – minimum voltage
Vmax (float) – maximum voltage
scan_speed (float) – Voltage scan speed [V/s]
direction (integer) – Perform a Vmin-Vmax-Vmin (1) or Vmax-Vmin-Vmax scan (-1)
steps (integer) – Number of time steps
G_frac (float) – Device Parameter | Fractional Generation rate
Vminexpo (float) – Voltage at which the exponential voltage steps start if Vmin or Vmax is 0
- Returns:
np.array – Array of time points
np.array – Array of voltages
np.array – Array of generation rates
- pySIMsalabim.experiments.hysteresis.calc_hysteresis_index(session_path, tj_file_name='tj.dat', tVG_file_ame='tVG.txt', plot_hyst_index=False)[source]
Calculate the hysteresis index from the simulated JV curve using the difference area between the forward and backward scan
- Parameters:
session_path (string) – working directory for zimt
tj_file_name (string) – Name of the tj file
tVG_file_ame (string) – Name of the tVG file
plot_hyst_index (bool) – If True, plot the JV curves with the normalisation and difference area
- Returns:
Hysteresis index
- Return type:
float
- pySIMsalabim.experiments.hysteresis.concatJVs(session_path, expJV_Vmin_Vmax, expJV_Vmax_Vmin, direction)[source]
Put the experimental forward and backward JV arrays together
- Parameters:
session_path (string) – working directory for zimt
expJV_Vmin_Vmax (string) – Name of the file of the forward JV scan
expJV_Vmax_Vmin (string) – Name of the file of the backward JV scan
direction (integer) – Perform a Vmin-Vmax-Vmin (1) or Vmax-Vmin-Vmax scan (-1)
- Returns:
Array of current and voltage of experimental JV
- Return type:
np.array
- pySIMsalabim.experiments.hysteresis.create_tVG_hysteresis(session_path, Vmin, Vmax, scan_speed, direction, steps, G_frac, tVG_name, **kwargs)[source]
Create a tVG file for hysteresis experiments.
- Parameters:
session_path (string) – working directory for zimt
Vmin (float) – Left voltage boundary
Vmax (float) – Right voltage boundary
scan_speed (float) – Voltage scan speed [V/s]
direction (integer) – Perform a Vmin-Vmax-Vmin (1) or Vmax-Vmin-Vmax scan (-1)
steps (integer) – Number of time steps
G_frac (float) – Device Parameter | Fractional Generation rate
tVG_name (string) – Device Parameter | Name of the tVG file
**kwargs (dict) – Additional keyword arguments
- Returns:
integer – Value to indicate the result of the process
string – A message to indicate the result of the process
- pySIMsalabim.experiments.hysteresis.plot_hysteresis_JV(session_path, path2file='tj.dat')[source]
Plot the hysteresis JV curve
- Parameters:
session_path (string) – working directory for zimt
path2file (string) – Path to the tj file
- Returns:
Axes object for the plot
- Return type:
Axes
- pySIMsalabim.experiments.hysteresis.read_Exp_JV(session_path, expJV_Vmin_Vmax, expJV_Vmax_Vmin)[source]
Read experimental forward and backward JV files
- Parameters:
session_path (string) – working directory for zimt
expJV_Vmin_Vmax (string) – Name of the file of the forward JV scan
expJV_Vmax_Vmin (string) – Name of the file of the backward JV scan
- Returns:
np.array – Array of current and voltage of experimental JV from Vmin to Vmax
np.array – Array of current and voltage of experimental JV from Vmax to Vmin
- pySIMsalabim.experiments.hysteresis.read_tj_file(session_path, tj_file_name='tj.dat')[source]
Read relevant parameters for admittance of the tj file
- Parameters:
session_path (string) – Path of the simulation folder for this session
data_tj (dataFrame) – Pandas dataFrame containing the tj output file from ZimT
- Returns:
Pandas dataFrame of the tj_file containing the time, current density, numerical error in the current density and the photogenerated current density
- Return type:
DataFrame
- pySIMsalabim.experiments.hysteresis.tVG_exp(session_path, expJV_Vmin_Vmax, expJV_Vmax_Vmin, scan_speed, direction, G_frac, tVG_name)[source]
Create a tVG file for hysteresis experiments where Vext is the same as the voltages in the experimental JV file
- Parameters:
session_path (string) – working directory for zimt
expJV_Vmin_Vmax (string) – Name of the file of the Vmin-Vmax JV scan
expJV_Vmax_Vmin (string) – Name of the file of the Vmax-Vmin JV scan
scan_speed (float) – Voltage scan speed [V/s]
direction (integer) – Perform a Vmin-Vmax-Vmin (1) or Vmax-Vmin-Vmax scan (-1)
G_frac (float) – Fractional Generation rate
tVG_name (string) – Name of the tVG file
- Returns:
integer – Value to indicate the result of the process
string – A message to indicate the result of the process
pySIMsalabim.experiments.impedance module
Perform impedance simulations
- pySIMsalabim.experiments.impedance.Bode_plot(session_path, output_file, xscale='log', yscale_1='linear', yscale_2='linear', plot_type=<function errorbar>)[source]
Plot the real and imaginary part of the impedance against frequency
- Parameters:
session_path (string) – working directory for zimt
xscale (string) – Scale of the x-axis. E.g linear or log
yscale_ax1 (string) – Scale of the left y-axis. E.g linear or log
yscale_ax2 (string) – Scale of the right y-axis. E.g linear or log
- pySIMsalabim.experiments.impedance.Capacitance_plot(session_path, output_file, xscale='log', yscale='linear')[source]
Plot the capacitance against frequency
- Parameters:
session_path (string) – working directory for zimt
xscale (string) – Scale of the x-axis. E.g linear or log
yscale (string) – Scale of the y-axis. E.g linear or log
- pySIMsalabim.experiments.impedance.Nyquist_plot(session_path, output_file, xscale='linear', yscale='linear', plot_type=<function errorbar>)[source]
Plot the real and imaginary part of the impedance against each other
- Parameters:
session_path (string) – working directory for zimt
xscale (string) – Scale of the x-axis. E.g linear or log
yscale (string) – Scale of the y-axis. E.g linear or log
- pySIMsalabim.experiments.impedance.calc_impedance(data, del_V, isToPlot, session_path, zimt_device_parameters, Rseries=0, Rshunt=-1000.0)[source]
Calculate the impedance over the frequency range
- Parameters:
data (dataFrame) – Pandas dataFrame containing the time, voltage, current density and numerical error in the current density of the tj_file
del_V (float) – Voltage step
isToPlot (list) – List of array indices that will be used in the plotting
- Returns:
np.array – Array of frequencies
np.array – Array of the real component of impedance
np.array – Array of the imaginary component of impedance
np.array – Array of complex error
np.array – Array of capacitance
np.array – Array of conductance
np.array – Array of error in capacitance
np.array – Array of error in conductance
- pySIMsalabim.experiments.impedance.calc_impedance_limit_time(I, errI, time, VStep, imax)[source]
Fourier Decomposition formula which computes the impedance at frequency freq (Hz) and its complex error Based on S.E. Laux, IEEE Trans. Electron Dev. 32 (10), 2028 (1985), eq. 5a, 5b We integrate from 0 to time[imax] at frequency 1/time[imax]
- Parameters:
I (np.array) – Array of currents
errI (np.array) – Numerical error in calculated currents (output of ZimT)
time (np.array) – Array with all time positions, from 0 to tmax
VStep (float) – Voltage step
imax (integer) – Index of the last timestep/first frequency for which the integrals are calculated
- Returns:
float – Frequency belonging to Z(f)
complex number – Impedance at frequency f: Z(f)
complex number – Numerical error in calculated impedance Z(f)
- pySIMsalabim.experiments.impedance.create_tVG_SS(V_0, G_frac, tVG_name, session_path)[source]
Creates the tVG file for the steady state simulation with only t 0 and V_0
- Parameters:
V_0 (float) – Voltage at t=0
del_V (float) – Voltage step that is applied after t=0
G_frac (float) – Fractional light intensity
tVG_name (string) – Name of the tVG file
session_path (string) – Path of the simulation folder for this session
- Returns:
A message to indicate the result of the process
- Return type:
string
- pySIMsalabim.experiments.impedance.create_tVG_impedance(V_0, del_V, G_frac, tVG_name, session_path, f_min, f_max, ini_timeFactor, timeFactor)[source]
Create a tVG file for impedance experiments.
- Parameters:
V_0 (float) – Voltage at t=0
del_V (float) – Voltage step that is applied after t=0
G_frac (float) – Fractional light intensity
tVG_name (string) – Name of the tVG file
session_path (string) – Path of the simulation folder for this session
f_min (float) – Minimum frequency
f_max (float) – Maximum frequency
ini_timeFactor (float) – Constant defining the size of the initial timestep
timeFactor (float) – Exponential increase of the timestep, to reduce the amount of timepoints necessary. Use values close to 1.
- Returns:
A message to indicate the result of the process
- Return type:
string
- pySIMsalabim.experiments.impedance.create_tVG_tolDens(V_0, del_V, ini_timeFactor, f_min, f_max, G, tVG_name, session_path)[source]
Creates the tVG file for the steady state simulation with only t 0 and V_0
- Parameters:
V_0 (float) – Voltage at t=0
del_V (float) – Voltage step that is applied after t=0
G_frac (float) – Fractional light intensity
tVG_name (string) – Name of the tVG file
session_path (string) – Path of the simulation folder for this session
- Returns:
A message to indicate the result of the process
- Return type:
string
- pySIMsalabim.experiments.impedance.get_impedance(data, f_min, f_max, f_steps, del_V, session_path, output_file, zimt_device_parameters, Rseries=0, Rshunt=-1000.0)[source]
Calculate the impedance from the simulation result
- Parameters:
data (DataFrame) – DataFrame with the simulation results (tj.dat) file
f_min (float) – Minimum frequency
f_max (float) – Maximum frequency
f_steps (float) – Frequency step
del_V (float) – Voltage step
session_path (string) – working directory for zimt
output_file (string) – name of the file where the impedance data is stored
- Returns:
returns -1 (failed) or 1 (success), including a message
- Return type:
integer,string
- pySIMsalabim.experiments.impedance.get_tolDens(zimt_device_parameters, session_path, f_min, f_max, V_0, G_frac, del_V, run_mode, tVG_name, tj_name, varFile, ini_timeFactor, dum_str, cmd_pars)[source]
Calculate the tolerance of the density solver, to ensure a reliable impedance spectrum
- Parameters:
zimt_device_parameters (string) – Name of the zimt device parameters file
session_path (string) – Working directory for zimt
f_min (float) – Minimum frequency
f_max (float) – Maximum frequency
V_0 (float) – Voltage at t=0
G_frac (float) – Fractional light intensity
del_V (float) – Voltage step
run_mode (bool) – Indicate whether the script is in ‘web’ mode (True) or standalone mode (False). Used to control the console output
tVG_name (string) – Name of the tVG file
tj_name (string) – Name of the tj file
varFile (string) – Name of the var file
ini_timeFactor (float) – Constant defining the size of the initial timestep
dum_str (string) – dummy string with UUID string to append to the file names
cmd_pars (list) – List of dictionaries with the command line parameters
- Returns:
integer – Return code of the simulation, 0 if successful, -1 if one of the steps failed else the return code of the simulation
string – Return message to display on the UI in case of failure
float – Tolerance of the density solver. If failed, returns None
- pySIMsalabim.experiments.impedance.plot_impedance(session_path, output_file='freqZ.dat')[source]
Make a Bode and Nyquist plot of the impedance
- Parameters:
session_path (string) – working directory for zimt
- pySIMsalabim.experiments.impedance.run_impedance_simu(zimt_device_parameters, session_path, f_min, f_max, f_steps, V_0, G_frac=1, del_V=0.01, run_mode=False, tVG_name='tVG.txt', output_file='freqZ.dat', tj_name='tj.dat', varFile='none', ini_timeFactor=0.001, timeFactor=1.02, **kwargs)[source]
Create a tVG file and run ZimT with impedance device parameters
- Parameters:
zimt_device_parameters (string) – Name of the zimt device parameters file
session_path (string) – Working directory for zimt
f_min (float) – Minimum frequency
f_max (float) – Maximum frequency
f_steps (float) – Frequency step
V_0 (float) – Voltage at t=0
G_frac (float, optional) – Fractional light intensity, by default 1
del_V (float, optional) – Voltage step, by default 0.01
run_mode (bool, optional) – Indicate whether the script is in ‘web’ mode (True) or standalone mode (False). Used to control the console output, by default False
tVG_name (string, optional) – Name of the tVG file, by default tVG.txt
output_file (string, optional) – Name of the file where the impedance data is stored, by default freqZ.dat
tj_name (string, optional) – Name of the tj file where the impedance data is stored, by default tj.dat
varFile (string, optional) – Name of the var file, by default ‘none’
ini_timeFactor (float, optional) – Constant defining the size of the initial timestep, by default 1e-3
timeFactor (float, optional) – Exponential increase of the timestep, to reduce the amount of timepoints necessary. Use values close to 1., by default 1.02
- Returns:
CompletedProcess – Output object of with returncode and console output of the simulation
string – Return message to display on the UI, for both success and failed
- pySIMsalabim.experiments.impedance.store_impedance_data(session_path, freq, ReZ, ImZ, errZ, C, G, errC, errG, output_file)[source]
Save the frequency, real & imaginary part of the impedance & impedance error in one file called freqZ.dat
- Parameters:
session_path (string) – working directory for zimt
freq (np.array) – Array of frequencies
ReZ (np.array) – Array of the real component of impedance
ImZ (np.array) – Array of the imaginary component of impedance
errZ (np.array) – Array of complex error
C (np.array) – Array of capacitance
G (np.array) – Array of conductance
errC (np.array) – Array of error in capacitance
errG (np.array) – Array of error in conductance
pySIMsalabim.experiments.imps module
Perform Intensity-Modulated Photocurrent Spectroscopy (IMPS) simulations
- pySIMsalabim.experiments.imps.ColeCole_plot(session_path, output_file, xscale='linear', yscale='linear', plot_type=<function errorbar>)[source]
Plot the Cole-Cole plot with the real and imaginary part of the admittance against frequency
- Parameters:
session_path (string) – working directory for zimt
output_file (string) – Filename where the admittance data is stored
xscale (string) – Scale of the x-axis. E.g linear or log
yscale (string) – Scale of the y-axis. E.g linear or log
plot_type (matplotlib.pyplot) – Type of plot to display
- pySIMsalabim.experiments.imps.IMPS_plot(session_path, output_file, xscale='log', yscale='log', plot_type=<function plot>)[source]
Plot the real and imaginary part of the admittance against frequency
- Parameters:
session_path (string) – working directory for zimt
xscale (string) – Scale of the x-axis. E.g linear or log
yscale_ax1 (string) – Scale of the left y-axis. E.g linear or log
yscale_ax2 (string) – Scale of the right y-axis. E.g linear or log
- pySIMsalabim.experiments.imps.calc_IMPS(data, isToPlot)[source]
Calculate the admittance over the frequency range
- Parameters:
data (dataFrame) – Pandas dataFrame containing the time, current density, numerical error in the current density and the photogenerated current density of the tj_file
isToPlot (list) – List of array indices that will be used in the plotting
- Returns:
np.array – Array of frequencies
np.array – Array of the real component of admittance
np.array – Array of the imaginary component of admittance
np.array – Array of complex error
- pySIMsalabim.experiments.imps.calc_IMPS_limit_time(I, errI, time, imax)[source]
Fourier Decomposition formula which computes the admittance at frequency freq (Hz) and its complex error Based on S.E. Laux, IEEE Trans. Electron Dev. 32 (10), 2028 (1985), eq. 5a, 5b We integrate from 0 to time[imax] at frequency 1/time[imax]
- Parameters:
I (np.array) – Array of currents
errI (np.array) – Numerical error in calculated currents (output of ZimT)
time (np.array) – Array with all time positions, from 0 to tmax
imax (integer) – Index of the last timestep/first frequency for which the integrals are calculated
- Returns:
float – Frequency belonging to Y(f)
complex number – Admittance at frequency f: Y(f)
complex number – Numerical error in calculated admittance Y(f)
- pySIMsalabim.experiments.imps.create_tVG_IMPS(V, G_frac, del_G, tVG_name, session_path, f_min, f_max, ini_timeFactor, timeFactor)[source]
Create a tVG file for IMPS experiments.
- Parameters:
V (float) – Voltage, the voltage is constant over the whole time range
G_frac (float) – Fractional light intensity
del_G (float) – Applied generation rate increase at t=0
tVG_name (string) – Name of the tVG file
session_path (string) – Path of the simulation folder for this session
f_min (float) – Minimum frequency
f_max (float) – Maximum frequency
ini_timeFactor (float) – Constant defining the size of the initial timestep
timeFactor (float) – Exponential increase of the timestep, to reduce the amount of timepoints necessary. Use values close to 1.
- Returns:
A message to indicate the result of the process
- Return type:
string
- pySIMsalabim.experiments.imps.get_IMPS(data, f_min, f_max, f_steps, session_path, output_file)[source]
Calculate the IMPS from the simulation result
- Parameters:
data (DataFrame) – DataFrame with the simulation results (tj.dat) file
f_min (float) – Minimum frequency
f_max (float) – Maximum frequency
f_steps (float) – Frequency step
session_path (string) – working directory for zimt
output_file (string) – name of the file where the IMPS data is stored
- Returns:
returns -1 (failed) or 1 (success), including a message
- Return type:
integer,string
- pySIMsalabim.experiments.imps.plot_IMPS(session_path, output_file='freqY.dat')[source]
Plot admittance of IMPS
- Parameters:
session_path (string) – working directory for zimt
- pySIMsalabim.experiments.imps.run_IMPS_simu(zimt_device_parameters, session_path, f_min, f_max, f_steps, V, G_frac, GStep=0.05, run_mode=False, tVG_name='tVG.txt', output_file='freqY.dat', tj_name='tj.dat', varFile='none', ini_timeFactor=0.001, timeFactor=1.02, **kwargs)[source]
Create a tVG file and run ZimT with admittance device parameters
- Parameters:
zimt_device_parameters (string) – Name of the zimt device parameters file
session_path (string) – Path of the simulation folder for this session
f_min (float) – Minimum frequency
f_max (float) – Maximum frequency
f_steps (float) – Frequency step
V (float) – Voltage, the voltage is constant over the whole time range
G_frac (float) – Fractional light intensity
GStep (float, optional) – Applied generation rate increase at t=0, by default 0.05
run_mode (bool, optional) – Indicate whether the script is in ‘web’ mode (True) or standalone mode (False). Used to control the console output, by default False
tVG_name (string, optional) – Name of the tVG file, by default tVG.txt
output_file (string, optional) – Name of the file where the admittance data is stored, by default freqY.dat
tj_name (string, optional) – Name of the tj file where the admittance data is stored, by default tj.dat
varFile (string, optional) – Name of the var file, by default ‘none’
ini_timeFactor (float, optional) – Constant defining the size of the initial timestep, by default 1e-3
timeFactor (float, optional) – Exponential increase of the timestep, to reduce the amount of timepoints necessary. Use values close to 1., by default 1.02
- Returns:
CompletedProcess – Output object of with returncode and console output of the simulation
string – Return message to display on the UI, for both success and failed
- pySIMsalabim.experiments.imps.store_IMPS_data(session_path, freq, ReY, ImY, errY, output_file)[source]
Save the frequency, real & imaginary part of the admittance & its error in one file called freqY.dat
- Parameters:
session_path (string) – working directory for zimt
freq (np.array) – Array of frequencies
ReY (np.array) – Array of the real component of admittance
ImY (np.array) – Array of the imaginary component of admittance
errY (np.array) – Array of complex error in admittance