Coverage for pySDC/implementations/problem_classes/BuckConverter.py: 100%

47 statements  

« prev     ^ index     » next       coverage.py v7.6.7, created at 2024-11-16 14:51 +0000

1import numpy as np 

2 

3from pySDC.core.problem import Problem 

4from pySDC.implementations.datatype_classes.mesh import mesh, imex_mesh 

5 

6 

7class buck_converter(Problem): 

8 r""" 

9 Example implementing the model of a buck converter, which is also called a step-down converter. The converter has two different 

10 states and each of this state can be expressed as a nonhomogeneous linear system of ordinary differential equations (ODEs) 

11 

12 .. math:: 

13 \frac{d u(t)}{dt} = A_k u(t) + f_k (t) 

14 

15 for :math:`k=1,2`. The two states are the following. Define :math:`T_{sw}:=\frac{1}{f_{sw}}` as the switching period with 

16 switching frequency :math:`f_{sw}`. The duty cycle :math:`d` defines the period of how long the switches are in one state 

17 until they switch to the other state. Roughly saying, the duty cycle can be seen as a percentage [1]_. A duty cycle of one means 

18 that the switches are always in only one state. If :math:`0 \leq \frac{t}{T_{sw}} \bmod 1 \leq d` [2]_: 

19 

20 .. math:: 

21 \frac{d v_{C_1} (t)}{dt} = -\frac{1}{R_s C_1}v_{C_1} (t) - \frac{1}{C_1} i_{L_1} (t) + \frac{V_s}{R_s C_1}, 

22 

23 .. math:: 

24 \frac{d v_{C_2} (t)}{dt} = -\frac{1}{R_\ell C_2}v_{C_2} (t) + \frac{1}{C_2} i_{L_1} (t), 

25 

26 .. math:: 

27 \frac{d i_{L_1} (t)}{dt} = \frac{1}{L_1} v_{C_1} (t) - \frac{1}{L_1} v_{C_2} (t) - \frac{R_\pi}{L_1} i_{L_1} (t), 

28 

29 Otherwise, the equations are 

30 

31 .. math:: 

32 \frac{d v_{C_1} (t)}{dt} = -\frac{1}{R_s C_1}v_{C_1} (t) + \frac{V_s}{R_s C_1}, 

33 

34 .. math:: 

35 \frac{d v_{C_2} (t)}{dt} = -\frac{1}{R_\ell C_2}v_{C_2} (t) + \frac{1}{C_2} i_{L_1} (t), 

36 

37 .. math:: 

38 \frac{d i_{L_1} (t)}{dt} = \frac{R_\pi}{R_s L_1} v_{C_1} (t) - \frac{1}{L_1} v_{C_2} (t) - \frac{R_\pi V_s}{L_1 R_s}. 

39 

40 using an initial condition. 

41 

42 Parameters 

43 ---------- 

44 duty : float, optional 

45 Duty cycle :math:`d` between zero and one indicates the time period how long the converter stays on one switching 

46 state until it switches to the other state. 

47 fsw : int, optional 

48 Switching frequency, it is used to determine the number of time steps after the switching state is changed. 

49 Vs : float, optional 

50 Voltage at the voltage source :math:`V_s`. 

51 Rs : float, optional 

52 Resistance of the resistor :math:`R_s` at the voltage source. 

53 C1 : float, optional 

54 Capacitance of the capacitor :math:`C_1`. 

55 Rp : float, optional 

56 Resistance of the resistor in front of the inductor :math:`R_\pi`. 

57 L1 : float, optional 

58 Inductance of the inductor :math:`L_1`. 

59 C2 : float, optional 

60 Capacitance of the capacitor :math:`C_2`. 

61 Rl : float, optional 

62 Resistance of the resistor :math:`R_\pi` 

63 

64 Attributes 

65 ---------- 

66 A : np.2darray 

67 Coefficient matrix of the ODE system. 

68 

69 Note 

70 ---- 

71 The duty cycle needs to be a value in :math:`[0,1]`. 

72 

73 References 

74 ---------- 

75 .. [1] J. Sun. Pulse-Width Modulation. 25-61. Springer, (2012). 

76 .. [2] J. Gyselinck, C. Martis, R. V. Sabariego. Using dedicated time-domain basis functions for the simulation of 

77 pulse-width-modulation controlled devices - application to the steady-state regime of a buck converter. Electromotion (2013). 

78 """ 

79 

80 dtype_u = mesh 

81 dtype_f = imex_mesh 

82 

83 def __init__(self, duty=0.5, fsw=1e3, Vs=10.0, Rs=0.5, C1=1e-3, Rp=0.01, L1=1e-3, C2=1e-3, Rl=10): 

84 """Initialization routine""" 

85 

86 # invoke super init, passing number of dofs 

87 nvars = 3 

88 super().__init__(init=(nvars, None, np.dtype('float64'))) 

89 self._makeAttributeAndRegister( 

90 'nvars', 'duty', 'fsw', 'Vs', 'Rs', 'C1', 'Rp', 'L1', 'C2', 'Rl', localVars=locals(), readOnly=True 

91 ) 

92 

93 self.A = np.zeros((nvars, nvars)) 

94 

95 def eval_f(self, u, t): 

96 """ 

97 Routine to evaluate the right-hand side of the problem. 

98 

99 Parameters 

100 ---------- 

101 u : dtype_u 

102 Current values of the numerical solution. 

103 t : float 

104 Current time of the numerical solution is computed. 

105 

106 Returns 

107 ------- 

108 f : dtype_f 

109 The right-hand side of the problem. 

110 """ 

111 Tsw = 1 / self.fsw 

112 

113 f = self.dtype_f(self.init, val=0.0) 

114 f.impl[:] = self.A.dot(u) 

115 

116 if 0 <= ((t / Tsw) % 1) <= self.duty: 

117 f.expl[0] = self.Vs / (self.Rs * self.C1) 

118 f.expl[2] = 0 

119 

120 else: 

121 f.expl[0] = self.Vs / (self.Rs * self.C1) 

122 f.expl[2] = -(self.Rp * self.Vs) / (self.L1 * self.Rs) 

123 

124 return f 

125 

126 def solve_system(self, rhs, factor, u0, t): 

127 r""" 

128 Simple linear solver for :math:`(I-factor\cdot A)\vec{u}=\vec{rhs}`. 

129 

130 Parameters 

131 ---------- 

132 rhs : dtype_f 

133 Right-hand side for the linear system. 

134 factor : float 

135 Abbrev. for the local stepsize (or any other factor required). 

136 u0 : dtype_u 

137 Initial guess for the iterative solver. 

138 t : float 

139 Current time (e.g. for time-dependent BCs). 

140 

141 Returns 

142 ------- 

143 me : dtype_u 

144 The solution as mesh. 

145 """ 

146 Tsw = 1 / self.fsw 

147 self.A = np.zeros((3, 3)) 

148 

149 if 0 <= ((t / Tsw) % 1) <= self.duty: 

150 self.A[0, 0] = -1 / (self.C1 * self.Rs) 

151 self.A[0, 2] = -1 / self.C1 

152 

153 self.A[1, 1] = -1 / (self.C2 * self.Rl) 

154 self.A[1, 2] = 1 / self.C2 

155 

156 self.A[2, 0] = 1 / self.L1 

157 self.A[2, 1] = -1 / self.L1 

158 self.A[2, 2] = -self.Rp / self.L1 

159 

160 else: 

161 self.A[0, 0] = -1 / (self.C1 * self.Rs) 

162 

163 self.A[1, 1] = -1 / (self.C2 * self.Rl) 

164 self.A[1, 2] = 1 / self.C2 

165 

166 self.A[2, 0] = self.Rp / (self.L1 * self.Rs) 

167 self.A[2, 1] = -1 / self.L1 

168 

169 me = self.dtype_u(self.init) 

170 me[:] = np.linalg.solve(np.eye(self.nvars) - factor * self.A, rhs) 

171 return me 

172 

173 def u_exact(self, t): 

174 r""" 

175 Routine to compute the exact solution at time :math:`t`. 

176 

177 Parameters 

178 ---------- 

179 t : float 

180 Time of the exact solution. 

181 

182 Returns 

183 ------- 

184 me : dtype_u 

185 The exact solution. 

186 """ 

187 assert t == 0, 'ERROR: u_exact only valid for t=0' 

188 

189 me = self.dtype_u(self.init) 

190 

191 me[0] = 0.0 # v1 

192 me[1] = 0.0 # v2 

193 me[2] = 0.0 # p3 

194 

195 return me