1-phase VSI (PWM)

The figure below shows a full-bridge bi-polar PWM voltage source inverter supplying an $RL$-load with $R=10\,\Omega$, $L=10\,$mH, and $V_{dc}=400\,$V. The switches S1 and S2 are switched in a complementary fashion with sinusoidal pulse width modulation technique. The switch S1 and S4 are turned on simultaneously. Similarly, the switch S3 and S2 are turned on simultaneously. The modulating voltage $V_m(t)=0.8\,\sin (100\,\pi\,t)$, and the triangular carrier voltage, varying between $+1$ and $-1$, has a frequency of $2\,$kHz. Determine
  1. RMS value of the output voltage
  2. RMS values of the fundamental components of the load voltage and load current
  3. average current through the DC source
  4. RMS values of the two lowest-order harmonics of the DC source current (50 Hz, 100 Hz)
  5. average power delivered to the load
  6. dominant harmonics present in the load voltage and current
In [1]:
from IPython.display import Image
Image(filename =r'VSI_1ph_6_fig_1.png', width=550)
Out[1]:
No description has been provided for this image
In [2]:
# run this cell to view the circuit file.
%pycat VSI_1ph_6_orig.in

We now replace the strings such as \$Vdc, \$L, with the values of our choice by running the python script given below. It takes an existing circuit file VSI_1ph_6_orig.in and produces a new circuit file VSI_1ph_6.in, after replacing \$Vdc, \$L, etc. with values of our choice.

In [3]:
import gseim_calc as calc
s_Vdc = "400"
s_L = "10e-3"
s_R = "10"
s_f_carrier = "2e3"
s_M = "0.8"
s_dt_min = "1e-6"
s_dt_nrml = "10e-6"

f_sin = 50.0
s_f_sin = ("%11.4E"%(f_sin)).strip()

T = 1/f_sin
s_2T = ("%11.4E"%(2.0*T)).strip()

l = [
  ('$Vdc', s_Vdc),
  ('$L', s_L),
  ('$R', s_R),
  ('$f_carrier', s_f_carrier),
  ('$f_sin', s_f_sin),
  ('$2T', s_2T),
  ('$M', s_M),
  ('$dt_min', s_dt_min),
  ('$dt_nrml', s_dt_nrml)
]
calc.replace_strings_1("VSI_1ph_6_orig.in", "VSI_1ph_6.in", l)
print('VSI_1ph_6.in is ready for execution')
VSI_1ph_6.in is ready for execution
Execute the following cell to run GSEIM on VSI_1ph_6.in.
In [4]:
import os
import dos_unix
# uncomment for windows:
#dos_unix.d2u("VSI_1ph_6.in")
os.system('run_gseim VSI_1ph_6.in')
get_lib_elements: filename gseim_aux/xbe.aux
get_lib_elements: filename gseim_aux/ebe.aux
Circuit: filename = VSI_1ph_6.in
main: i_solve = 0
main: calling solve_trns
Transient simulation starts...
i=0
i=1000
i=2000
i=3000
i=4000
GSEIM: Program completed.
Out[4]:
0

The circuit file (VSI_1ph_6.in) is created in the same directory as that used for launching Jupyter notebook. The last step (i.e., running GSEIM on VSI_1ph_6.in) creates a data file called VSI_1ph_6.datin the same directory. We can now use the python code below to compute/plot the various quantities of interest.

In [5]:
import numpy as np
import matplotlib.pyplot as plt 
import gseim_calc as calc
from setsize import set_size

slv = calc.slv("VSI_1ph_6.in")

i_slv = 0
i_out = 0
filename = slv.l_filename_all[i_slv][i_out]
print('filename:', filename)
u = np.loadtxt(filename)
t = u[:, 0]

col_v_out = slv.get_index(i_slv,i_out,"v_out")
col_IR    = slv.get_index(i_slv,i_out,"IR"   )
col_ISrc  = slv.get_index(i_slv,i_out,"ISrc" )
col_P_R   = slv.get_index(i_slv,i_out,"P_R"  )
col_g1    = slv.get_index(i_slv,i_out,"g1"   )
col_g2    = slv.get_index(i_slv,i_out,"g2"   )
col_s     = slv.get_index(i_slv,i_out,"s"    )
col_t     = slv.get_index(i_slv,i_out,"t"    )

# since we have stored two cycles, we need to divide the last time point
# by 2 to get the period:

T = t[-1]/2

l_IR    = calc.avg_rms_2(t, u[:,col_IR],    T, 2.0*T, 1.0e-4*T)
l_v_out = calc.avg_rms_2(t, u[:,col_v_out], T, 2.0*T, 1.0e-4*T)
l_P_R   = calc.avg_rms_2(t, u[:,col_P_R],   T, 2.0*T, 1.0e-4*T)
l_ISrc  = calc.avg_rms_2(t, u[:,col_ISrc],  T, 2.0*T, 1.0e-4*T)

print('rms load voltage:', "%11.4E"%l_v_out[2][0])
print('average source current:', "%11.4E"%l_ISrc[1][0])
print('average power delivered to load:', "%11.4E"%l_P_R[1][0])

color1='green'
color2='crimson'
color3='goldenrod'
color4='blue'
color5='cornflowerblue'

fig, ax = plt.subplots(4, sharex=False)
plt.subplots_adjust(wspace=0, hspace=0.0)

set_size(7, 7, ax[0])

for i in range(4):
    ax[i].set_xlim(left=0, right=T*1e3)

ax[0].set_ylim(bottom=-1.4, top=1.4)
ax[1].set_ylim(bottom=-500, top=500)

ax[0].grid(color='#CCCCCC', linestyle='solid', linewidth=0.5)
ax[1].grid(color='#CCCCCC', linestyle='solid', linewidth=0.5)
ax[2].grid(color='#CCCCCC', linestyle='solid', linewidth=0.5)
ax[3].grid(color='#CCCCCC', linestyle='solid', linewidth=0.5)

ax[0].plot((t-T)*1e3, u[:,col_t    ], color=color5, linewidth=1.0, label="$t$")
ax[0].plot((t-T)*1e3, u[:,col_s    ], color=color4, linewidth=1.0, label="$s$")
ax[1].plot((t-T)*1e3, u[:,col_v_out], color=color1, linewidth=1.0, label="$V_{out}$")
ax[2].plot((t-T)*1e3, u[:,col_ISrc ], color=color3, linewidth=1.0, label="$i_{dc}$")
ax[3].plot((t-T)*1e3, u[:,col_IR   ], color=color2, linewidth=1.0, label="$i_L$")

ax[1].set_ylabel(r'$V_{out}$', fontsize=14)
ax[2].set_ylabel(r'$i_{dc}$',  fontsize=14)
ax[3].set_ylabel(r'$i_L$',     fontsize=14)

ax[3].set_xlabel('time (msec)', fontsize=14)

ax[0].legend(loc = 'lower right',frameon = True, fontsize = 10, title = None,
   markerfirst = True, markerscale = 1.0, labelspacing = 0.5, columnspacing = 2.0,
   prop = {'size' : 12},)

ax[0].tick_params(labelbottom=False)
ax[1].tick_params(labelbottom=False)
ax[2].tick_params(labelbottom=False)

#plt.tight_layout()
plt.show()
filename: VSI_1ph_6.dat
rms load voltage:  4.0000E+02
average source current:  1.1772E+01
average power delivered to load:  4.6969E+03
No description has been provided for this image
In [6]:
import numpy as np
import matplotlib.pyplot as plt 
from matplotlib.ticker import (MultipleLocator, AutoMinorLocator)
import gseim_calc as calc
from setsize import set_size

slv = calc.slv("VSI_1ph_6.in")

i_slv = 0
i_out = 0
filename = slv.l_filename_all[i_slv][i_out]
print('filename:', filename)
u = np.loadtxt(filename)
t = u[:, 0]

col_v_out = slv.get_index(i_slv,i_out,"v_out")
col_IR    = slv.get_index(i_slv,i_out,"IR"   )
col_ISrc  = slv.get_index(i_slv,i_out,"ISrc" )
col_P_R   = slv.get_index(i_slv,i_out,"P_R"  )
col_g1    = slv.get_index(i_slv,i_out,"g1"   )
col_g2    = slv.get_index(i_slv,i_out,"g2"   )
col_s     = slv.get_index(i_slv,i_out,"s"    )
col_t     = slv.get_index(i_slv,i_out,"t"    )

T = t[-1]/2

# compute Fourier coeffs:

t_start = T 
t_end = 2.0*T

n_fourier = 45

coeff_IR, thd_IR = calc.fourier_coeff_1C(t, u[:,col_IR], 
    t_start, t_end, 1.0e-8, n_fourier)

coeff_v_out, thd_v_out = calc.fourier_coeff_1C(t, u[:,col_v_out], 
    t_start, t_end, 1.0e-8, n_fourier)

coeff_ISrc, thd_ISrc = calc.fourier_coeff_1C(t, u[:,col_ISrc], 
    t_start, t_end, 1.0e-8, n_fourier)

print("THD (load current): ", "%5.2f"%(thd_IR*100.0), "%")
print("load current fundamental: RMS value: ", "%11.4E"%(coeff_IR[1]/np.sqrt(2.0)))
print("THD (load voltage): ", "%5.2f"%(thd_v_out*100.0), "%")
print("load voltage fundamental: RMS value: ", "%11.4E"%(coeff_v_out[1]/np.sqrt(2.0)))
print("DC source current fundamental: RMS value: ", "%11.4E"%(coeff_ISrc[1]/np.sqrt(2.0)))
print("DC source current 1st harmonic: RMS value: ", "%11.4E"%(coeff_ISrc[2]/np.sqrt(2.0)))

x = np.linspace(0, n_fourier, n_fourier+1)

y_IR    = np.array(coeff_IR)
y_v_out = np.array(coeff_v_out)
y_ISrc  = np.array(coeff_ISrc)

fig, ax = plt.subplots(3, sharex=False)
plt.subplots_adjust(wspace=0, hspace=0.0)
grid_color='#CCCCCC'

set_size(7, 6, ax[0])

delta = 5.0
x_major_ticks = np.arange(0.0, (float(n_fourier+1)), delta)
x_minor_ticks = np.arange(0.0, (float(n_fourier+1)), delta/5)

for i in range(3):
    ax[i].set_xlim(left=-1.0, right=float(n_fourier))
    ax[i].set_xticks(x_major_ticks)
    ax[i].set_xticks(x_minor_ticks, minor=True)
    ax[i].grid(visible=True, which='major', axis='x', color=grid_color, linestyle='-', zorder=0)
#   ax[i].grid(visible=True, which='minor', axis='x', color=grid_color, linestyle='-', zorder=0)

ax[0].set_ylabel('$i_{load}$',fontsize=14)
ax[1].set_ylabel('$v_{out}$', fontsize=14)
ax[2].set_ylabel('$i_{dc}$',  fontsize=14)

ax[2].set_xlabel('N', fontsize=14)

bars1 = ax[0].bar(x, y_IR,    width=0.3, color='red',   label="$i_{load}$", zorder=3)
bars2 = ax[1].bar(x, y_v_out, width=0.3, color='blue',  label="$V_{out}$",  zorder=3)
bars3 = ax[2].bar(x, y_ISrc,  width=0.3, color='green', label="$i_{dc}$",   zorder=3)

plt.tight_layout()
plt.show()
filename: VSI_1ph_6.dat
THD (load current):   9.54 %
load current fundamental: RMS value:   2.1574E+01
THD (load voltage):  145.78 %
load voltage fundamental: RMS value:   2.2627E+02
DC source current fundamental: RMS value:   2.9626E-03
DC source current 1st harmonic: RMS value:   8.5677E+00
No description has been provided for this image

This notebook was contributed by Prof. Nakul Narayanan K, Govt. Engineering College, Thrissur. He may be contacted at nakul@gectcr.ac.in.

In [ ]: