rf_control#

Design and analyze RF feedback/feedforward controllers.

Functions

ADRC_control_step(Akd, Bkd, Ckd, Dkd, Aobd, ...)

Controller execute for one step based on the controller's discrete state-space equation including the ADRC observer.

ADRC_controller(half_bw[, pole_scale])

Generate the continous ADRC controller.

AFF_ilc(vc_err, L)

Apply the ILC.

AFF_ilc_design(h, pulw[, P, Q])

Adaptive feedforward with optimal iterative learning control (ILC).

AFF_timerev_lpf(vfb, fcut, fs[, vff_cor])

Time-reversed low-pass filter, we only apply the first order IIR low-pass, which gives up to 90 degrees phase lead.

basic_rf_controller(Kp, Ki[, notch_conf, ...])

Generate the continous state-space equation for a basic controller with

cav_sp_ff(half_bw, filling_len, flattop_len, ...)

Generate the basic cavity set point and feedforward used to configure the LLRF controller (later may apply smooth in the edges).

control_step(Akd, Bkd, Ckd, Dkd, err_step, ...)

Controller execute for one step based on the discrete state-space equation.

loop_analysis(AG, BG, CG, DG, AK, BK, CK, DK)

Control loop analysis, including

resp_inv_lsm(R[, regu])

Response matrix inversion with least-square method.

resp_inv_svd(R[, singular_val_filt])

Response matrix inversion with SVD.

ss_cascade(A1, B1, C1, D1, A2, B2, C2, D2[, Ts])

Cascade two state-space systems with the first system applyed to the input first.

ss_discrete(Ac, Bc, Cc, Dc, Ts[, method, ...])

Derive the discrete state-space equation from a continous one and compare their frequency responses.

ss_freqresp(A, B, C, D[, Ts, plot, ...])

Plot the frequency response of a state-space system.

ADRC_control_step(Akd, Bkd, Ckd, Dkd, Aobd, Bobd, b0, sp_step, vc_step, vd_step, state_k0, state_ob0, vf_step=0.0, ff_step=0.0, apply_to_err=False)[source]#

Controller execute for one step based on the controller’s discrete state-space equation including the ADRC observer.

Refer to LLRF Book Figure 4.6.

Parameters:
  • Akd – numpy matrix (complex), discrete RF controller

  • Bkd – numpy matrix (complex), discrete RF controller

  • Ckd – numpy matrix (complex), discrete RF controller

  • Dkd – numpy matrix (complex), discrete RF controller

  • Aobd – numpy matrix (complex), discrete ADRC observer (C,D not needed)

  • Bobd – numpy matrix (complex), discrete ADRC observer (C,D not needed)

  • b0 – float, ADRC gain

  • sp_step – complex, cavity voltage setpoint of this time step

  • vc_step – complex, cavity voltage meas. of last time step

  • vd_step – complex, cavity drive of last time step

  • vf_step – complex, feedfoward part of cavity drive of last time step

  • state_k0 – numpy matrix (complex), controller state of last time step

  • state_ob0 – numpy matrix (complex), ADRC observer state of last time step

  • ff_step – complex, feedforward of this time step

  • apply_to_err – boolean, True to apply ADRC to error, or apply to whole cavity voltage

Returns:
  • status – boolean, success (True) or fail (False)

  • ctrl_step – complex, overall output of the controller of this step

  • ctrl_out – complex, feedback output (exclude feedforward wrt ctrl_step)

  • state_k – numpy matrix (complex), controller state of this step (should input to the next execution of this function)

  • state_ob – numpy matrix (complex), ADRC observer state of this step (should input to the next execution of this function)

  • f – complex, estimated general disturbance by the ADRC observer

Question:

In the second equation, shall we use state_k0 or state_k?

Note

  1. looks like ADRC can mitigate the instability caused by the cavity’s passband modes, then no notch filter is required to filter the cavity voltage measurement.

  2. basic ADRC is perfect with a proportional controller - resulting in zero steady-state error. Looks like with other controller (e.g., PI control), ADRC needs to be revised for good performance.

  3. adding feedforward to the basic ADRC also causes some problem, like the steady state error becomes nonzero.

  4. The solution to handle the problems of point 2 and 3 above is to apply ADRC to the errors of the cavity drive and voltage. See the paper: S. Zhao et al, Tracking and disturbance rejection in non-minimum phase systems, Proceedings of the 33rd Chinese Control Conference, pp. 3834-3839, July 28-30, 2014, Nanjing, China.

ADRC_controller(half_bw, pole_scale=50.0)[source]#

Generate the continous ADRC controller. We assume that the system gain has been normalized to 1 and the system phase corrected to 0. Therefore, the referred cavity equation is dvc/dt + (half_bw - 1j * detuning)*vc = half_bw * vd where vc is the vector of cavity voltage and vd is the vector of cavity drive (in principle vd = 2 * beta * vf / (beta + 1)).

Refer to LLRF Book section 4.2.3.

Parameters:
  • half_bw – float, half bandwidth of the cavity, rad/s

  • pole_scale – float, define the pole location of the observer

Returns:
  • status – boolean, success (True) or fail (False)

  • Aobc, Bobc, Cobc, Dobc – numpy matrix (complex), continous ADRC observer

  • b0 – float, ADRC gain

AFF_ilc(vc_err, L)[source]#

Apply the ILC.

Refer to LLRF Book section 4.5.2.

Parameters:
  • vc_err – numpy array (complex), error of the cavity voltage waveform

  • L – numpy matrix (complex), gain matrix of ILC

Returns:

vff_cor – numpy array (complex), feedforward correction waveform

AFF_ilc_design(h, pulw, P=None, Q=None)[source]#

Adaptive feedforward with optimal iterative learning control (ILC).

Refer to LLRF Book section 4.5.2.

Parameters:
  • h – numpy array (complex), impulse response

  • pulw – int, pulse width as number of points

  • P – numpy matrix, positive-definite weight matrices

  • Q – numpy matrix, positive-definite weight matrices

Returns:
  • status – boolean, success (True) or fail (False)

  • L – numpy matrix (complex), gain matrix of ILC

AFF_timerev_lpf(vfb, fcut, fs, vff_cor=None)[source]#

Time-reversed low-pass filter, we only apply the first order IIR low-pass, which gives up to 90 degrees phase lead.

Refer to LLRF Book section 4.5.1.

Parameters:
  • vfb – numpy array (complex), feedback control waveform

  • fcut – float, cut-off frequency of the low-pass filter, Hz

  • fs – float, sampling frequency, Hz

  • vff_cor – numpy array (complex), buffer storing the filtered waveform

Returns:
  • status – boolean, success (True) or fail (False)

  • vff_cor – numpy array (complex), feedforward correction waveform

basic_rf_controller(Kp, Ki, notch_conf=None, plot=False, plot_pno=1000, plot_maxf=0.0)[source]#
Generate the continous state-space equation for a basic controller with
  • I/Q control strategy (input/output are all complex signals).

  • configurable for PI control + frequency notches.

Parameters:
  • Kp – float, proportional (P) gain

  • Ki – float, integral (I) gain

  • notch_conf – dict, with the following items: freq_offs: list, offset frequency to be notched, Hz; gain: list, gain of the notch (inverse of suppress ratio); half_bw: list, half bandwidth of the notch filter, rad/s

  • plot – boolean, enable the plot of frequency response of the controller

  • plot_pno – int, number of point in the plot

  • plot_maxf – float, frequency range (+-) to be plotted, Hz

Returns:
  • status – boolean, success (True) or fail (False)

  • Akc, Bkc, Ckc, Dkc – numpy matrix (complex), continuous state-space controller

cav_sp_ff(half_bw, filling_len, flattop_len, Ts, pno, vc0=1.0, detuning=0.0, beta=10000.0, const_fpow=True, ib0=None, phib_deg=0.0, beam_ids=0, beam_ide=0, roQ_or_RoQ=0.0, QL=3000000.0, machine='linac')[source]#

Generate the basic cavity set point and feedforward used to configure the LLRF controller (later may apply smooth in the edges).

Refer to LLRF Book section 9.2.1.

Parameters:
  • half_bw – float, half bandwidth of the cavity, rad/s

  • filling_len – int, length of cavity filling time, number of samples

  • flattop_len – int, length of the flattop for beam acc., number of samples

  • Ts – float, sampling time, s

  • pno – int, number of samples in the returned waveforms

  • vc0 – float, desired cavity voltage at the flattop, V

  • detuning – float, detuning of the cavity, rad/s

  • beta – float, input coupling factor (needed for NC cavities; for SC cavities, can use the default value, or you can specify it if more accurate result is needed)

  • const_fpow – boolean, True for forcing constant filling drive power/phase

  • ib0 – float, average beam current, A

  • phib_deg – float, beam accelerating phase (0 for on-crest), deg

  • beam_ids – int, beam starting sample index

  • beam_ide – int, beam ending sample index

  • roQ_or_RoQ – r/Q of Linac or R/Q of circular accelerator, Ohm (see the note below)

  • QL – float, loaded quality factor of the cavity

  • machine – string, linac or circular, used to select r/Q or R/Q

Returns:
  • status – boolean, success (True) or fail (False)

  • sp – numpy array (complex), set point waveform (for controller)

  • vf_ff – numpy array (complex), feedforward waveform (for controller)

  • vb – numpy array (complex), beam drive voltage waveform

  • Tpul – numpy array, time array for the waveforms

Note

Linacs define the r/Q = Vacc**2 / (w0 * U) while circular machines define R/Q = Vacc**2 / (2 * w0 * U), where Vacc is the accelerating voltage, w0 is the angular cavity resonance frequency and U is the cavity energy storage. Therefore, to use this function, one needs to specify the machine to be linac or circular. Generally, we have R/Q = 1/2 * r/Q.

control_step(Akd, Bkd, Ckd, Dkd, err_step, state_k0, ff_step=0.0)[source]#

Controller execute for one step based on the discrete state-space equation.

Parameters:
  • Akd – numpy matrix (complex), discrete state-space controller

  • Bkd – numpy matrix (complex), discrete state-space controller

  • Ckd – numpy matrix (complex), discrete state-space controller

  • Dkd – numpy matrix (complex), discrete state-space controller

  • err_step – complex, system output error (input to controller) of this step

  • state_k0 – numpy matrix (complex), state of the last step

  • ff_step – complex, feedforward of this step

Returns:
  • status – boolean, success (True) or fail (False)

  • ctrl_step – complex, overall output of the controller of this step

  • ctrl_out – complex, feedback output (exclude feedforward wrt ctrl_step)

  • state_k – numpy matrix (complex), state of this step (should be input to the next execution of this function)

Question:

In the second equation, shall we use state_k0 or state_k?

loop_analysis(AG, BG, CG, DG, AK, BK, CK, DK, Ts=None, delay_s=0, plot=True, plot_pno=100000, plot_maxf=0.0, label='')[source]#
Control loop analysis, including
  • derive the open loop transfer function.

  • calculate the senstivity and complementary sensitivity functions.

This function works for both continous system (Ts is None) and discrete systems (Ts has a nonzero floating value).

Parameters:
  • AG – numpy matrix (complex), plant model

  • BG – numpy matrix (complex), plant model

  • CG – numpy matrix (complex), plant model

  • DG – numpy matrix (complex), plant model

  • AK – numpy matrix (complex), controller

  • BK – numpy matrix (complex), controller

  • CK – numpy matrix (complex), controller

  • DK – numpy matrix (complex), controller

  • Ts – float, sampling time, s

  • delay_s – float, loop delay, s

  • plot – boolean, enable the plot of bode and Nyquist plots

  • plot_pno – int, number of point in the plot

  • plot_maxf – float, frequency range (+-) to be plotted, Hz

Returns:
  • status – boolean, success (True) or fail (False)

  • S_max – float, maximum sensitivity, dB

  • T_max – float, maximum complementary sensitivity, dB

resp_inv_lsm(R, regu=0.0)[source]#

Response matrix inversion with least-square method.

Refer to Beam Control Book section 2.4.3.

Parameters:
  • R – numpy matrix, response matrix

  • regu – float, regularization factor (should > 0)

Returns:

Rinv – numpy matrix, inversion of the response matrix

resp_inv_svd(R, singular_val_filt=0.0)[source]#

Response matrix inversion with SVD.

Refer to Beam Control Book section 2.4.2.

Parameters:
  • R – numpy matrix, response matrix

  • singular_val_filt – float, threshold of singular values, the ones smaller or equal to it will be discarded

Returns:

Rinv – numpy matrix, inversion of the response matrix

ss_cascade(A1, B1, C1, D1, A2, B2, C2, D2, Ts=None)[source]#

Cascade two state-space systems with the first system applyed to the input first. This function works for both continous system (Ts is None) and discrete systems (Ts has a nonzero floating value).

Parameters:
  • A1 – numpy matrix (complex), state-space model of system 1

  • B1 – numpy matrix (complex), state-space model of system 1

  • C1 – numpy matrix (complex), state-space model of system 1

  • D1 – numpy matrix (complex), state-space model of system 1

  • A2 – numpy matrix (complex), state-space model of system 2

  • B2 – numpy matrix (complex), state-space model of system 2

  • C2 – numpy matrix (complex), state-space model of system 2

  • D2 – numpy matrix (complex), state-space model of system 2

  • Ts – float, sampling time, s

Returns:
  • status boolean, success (True) or fail (False)

  • A, B, C, D – numpy matrix (complex), state-space model of cascaded system

ss_discrete(Ac, Bc, Cc, Dc, Ts, method='zoh', alpha=0.3, plot=False, plot_pno=1000, spec_data=False)[source]#

Derive the discrete state-space equation from a continous one and compare their frequency responses.

Parameters:
  • Ac – numpy matrix (complex), continous state-space model

  • Bc – numpy matrix (complex), continous state-space model

  • Cc – numpy matrix (complex), continous state-space model

  • Dc – numpy matrix (complex), continous state-space model

  • Ts – float, sampling time, s

  • method – string, gbt, bilinear, euler, backward_diff, zoh, foh or impulse (see document of signal.cont2discrete)

  • alpha – float, 0 to 1 (see document of signal.cont2discrete)

  • plot – boolean, enable the plot of frequency responses

  • plot_pno – int, number of point in the plot

  • spec_data – boolean, return spectra data or not

Returns:
  • status – boolean, success (True) or fail (False)

  • Ad, Bd, Cd, Dd – numpy matrix (complex), discrete state-space model

ss_freqresp(A, B, C, D, Ts=None, plot=False, plot_pno=1000, plot_maxf=0.0, title='Frequency Response')[source]#

Plot the frequency response of a state-space system. This function works for both continous system (Ts is None) and discrete systems (Ts has a nonzero floating value).

Parameters:
  • A – numpy matrix (complex), state-space model of system

  • B – numpy matrix (complex), state-space model of system

  • C – numpy matrix (complex), state-space model of system

  • D – numpy matrix (complex), state-space model of system

  • Ts – float, sampling time, s

  • plot – boolean, enable the plot of frequency response

  • plot_pno – int, number of point in the plot

  • plot_maxf – float, frequency range (+-) to be plotted, Hz

  • title – string, title showed on the plot

Returns:
  • status – boolean, success (True) or fail (False)

  • f_wf – numpy array, frequency waveform, Hz

  • A_wf_dB – numpy array, amplitude response waveform, dB

  • P_wf_deg – numpy array, phase response waveform, deg

  • h – numpy array (complex), complex response