implementations.problem_classes.Piline module

class piline(Vs=100.0, Rs=1.0, C1=1.0, Rpi=0.2, Lpi=1.0, C2=1.0, Rl=5.0)[source]

Bases: ptype

Example implementing the model of the piline. It serves as a transmission line in an energy grid. The problem of simulating the piline consists of three ordinary differential equations (ODEs) with nonhomogeneous part:

\[\frac{d v_{C_1} (t)}{dt} = -\frac{1}{R_s C_1}v_{C_1} (t) - \frac{1}{C_1} i_{L_\pi} (t) + \frac{V_s}{R_s C_1},\]
\[\frac{d v_{C_2} (t)}{dt} = -\frac{1}{R_\ell C_2}v_{C_2} (t) + \frac{1}{C_2} i_{L_\pi} (t),\]
\[\frac{d i_{L_\pi} (t)}{dt} = \frac{1}{L_\pi} v_{C_1} (t) - \frac{1}{L_\pi} v_{C_2} (t) - \frac{R_\pi}{L_\pi} i_{L_\pi} (t),\]

which can be expressed as a nonhomogeneous linear system of ODEs

\[\frac{d u(t)}{dt} = A u(t) + f(t)\]

using an initial condition.

Parameters:
  • Vs (float, optional) – Voltage at the voltage source \(V_s\).

  • Rs (float, optional) – Resistance of the resistor \(R_s\) at the voltage source.

  • C1 (float, optional) – Capacitance of the capacitor \(C_1\).

  • Rpi (float, optional) – Resistance of the resistor \(R_\pi\).

  • Lpi (float, optional) – Inductance of the inductor \(L_\pi\).

  • C2 (float, optional) – Capacitance of the capacitor \(C_2\).

  • Rl (float, optional) – Resistance of the resistive load \(R_\ell\).

  • Attributes

  • A (np.2darray) – Coefficient matrix of the linear ODE system.

dtype_f

alias of imex_mesh

dtype_u

alias of mesh

eval_f(u, t)[source]

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

Parameters:
  • u (dtype_u) – Current values of the numerical solution.

  • t (float) – Current time of the numerical solution is computed.

Returns:

f – The right-hand side of the problem.

Return type:

dtype_f

solve_system(rhs, factor, u0, t)[source]

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

Parameters:
  • rhs (dtype_f) – Right-hand side for the linear system.

  • factor (float) – Abbrev. for the local stepsize (or any other factor required).

  • u0 (dtype_u) – Initial guess for the iterative solver.

  • t (float) – Current time (e.g. for time-dependent BCs).

Returns:

me – The solution as mesh.

Return type:

dtype_u

u_exact(t, u_init=None, t_init=None)[source]

Routine to approximate the exact solution at time \(t\) by SciPy as a reference.

Parameters:
  • t (float) – Time of the exact solution.

  • u_init (pySDC.problem.Piline.dtype_u) – Initial conditions for getting the exact solution.

  • t_init (float) – The starting time.

Returns:

me – The reference solution.

Return type:

dtype_u