mbtrack2.tracking.rf module

This module handles radio-frequency (RF) cavitiy elements.

class RFCavity(ring: Synchrotron, m: int, Vc: float, theta: float)[source]

Bases: Element

Perfect RF cavity class for main and harmonic RF cavities. Use cosine definition.

Parameters

ring : Synchrotron object m : int

Harmonic number of the cavity

Vcfloat

Amplitude of cavity voltage [V]

thetafloat

Phase of Cavity voltage

__init__(ring: Synchrotron, m: int, Vc: float, theta: float)[source]
track(bunch: Bunch | Beam)[source]

Tracking method for the element. No bunch to bunch interaction, so written for Bunch objects and @Element.parallel is used to handle Beam objects.

Parameters

bunch : Bunch or Beam object

value(val: ndarray[tuple[int, ...], dtype[_ScalarType_co]]) ndarray[tuple[int, ...], dtype[_ScalarType_co]][source]
_abc_impl = <_abc._abc_data object>
class CavityResonator(ring: Synchrotron, m: int | float, Rs: float, Q: float, QL: float, detune: float, Ncav: int = 1, Vc: float = 0, theta: float = 0, n_bin: int = 75)[source]

Bases: object

Cavity resonator class for active or passive RF cavity with beam loading or HOM, based on [1,2].

Use cosine definition.

If used with mpi, beam.mpi.share_distributions must be called before the track method call.

Different kind of RF feeback and loops can be added using:

cavity_resonator.feedback.append(loop)

Parameters

ring : Synchrotron object m : int or float

Harmonic number of the cavity.

Rsfloat

Shunt impedance of the cavities in [Ohm], defined as 0.5*Vc*Vc/Pc. If Ncav = 1, used for the total shunt impedance. If Ncav > 1, used for the shunt impedance per cavity.

Qfloat

Quality factor of the cavity.

QLfloat

Loaded quality factor of the cavity.

detunefloat

Detuing of the cavity in [Hz], defined as (fr - m*ring.f1).

Ncavint, optional

Number of cavities.

Vcfloat, optinal

Total cavity voltage objective value in [V].

thetafloat, optional

Total cavity phase objective value in [rad].

n_binint, optional

Number of bins used for the beam loading computation. Only used if MPI is not used, otherwise n_bin must be specified in the beam.mpi.share_distributions method. The default is 75.

Attributes

beam_phasorcomplex

Beam phasor in [V].

beam_phasor_recordarray of complex

Last beam phasor value of each bunch in [V].

generator_phasorcomplex

Generator phasor in [V].

generator_phasor_recordarray of complex

Last generator phasor value of each bunch in [V].

cavity_phasorcomplex

Cavity phasor in [V].

cavity_phasor_recordarray of complex

Last cavity phasor value of each bunch in [V].

ig_phasor_recordarray of complex

Last current generator phasor of each bunch in [A]. Only used for some feedback types.

cavity_voltagefloat

Cavity total voltage in [V].

cavity_phasefloat

Cavity total phase in [rad].

loss_factorfloat

Cavity loss factor in [V/C].

Rs_per_cavityfloat

Shunt impedance of a single cavity in [Ohm], defined as 0.5*Vc*Vc/Pc.

betafloat

Coupling coefficient of the cavity.

frfloat

Resonance frequency of the cavity in [Hz].

wrfloat

Angular resonance frequency in [Hz.rad].

psifloat

Tuning angle in [rad].

filling_timefloat

Cavity filling time in [s].

Pcfloat

Power dissipated in the cavity walls in [W].

Pgfloat

Generator power in [W].

Vgrfloat

Generator voltage at resonance in [V].

theta_grfloat

Generator phase at resonance in [rad].

Vgfloat

Generator voltage in [V].

theta_gfloat

Generator phase in [rad].

trackingbool

True if the tracking has been initialized.

bunch_indexint

Number of the tracked bunch in the current core.

distancearray

Distance between bunches.

valid_bunch_index : array

Methods

Vbr(I0)

Return beam voltage at resonance in [V].

Vb(I0)

Return beam voltage in [V].

Pb(I0)

Return power transmitted to the beam in [W].

Pr(I0)

Return power reflected back to the generator in [W].

Z(f)

Cavity impedance in [Ohm] for a given frequency f in [Hz].

set_optimal_coupling(I0)

Set coupling to optimal value.

set_optimal_detune(I0)

Set detuning to optimal conditions.

set_generator(I0)

Set generator parameters.

plot_phasor(I0)

Plot phasor diagram.

is_DC_Robinson_stable(I0)

Check DC Robinson stability.

is_CBI_stable(I0)

Check Coupled-Bunch-Instability stability.

plot_DC_Robinson_stability()

Plot DC Robinson stability limit.

init_tracking(beam)

Initialization of the tracking.

track(beam)

Tracking method.

phasor_decay(time)

Compute the beam phasor decay during a given time span.

phasor_evol(profile, bin_length, charge_per_mp)

Compute the beam phasor evolution during the crossing of a bunch.

VRF(z, I0)

Return the total RF voltage.

dVRF(z, I0)

Return derivative of total RF voltage.

ddVRF(z, I0)

Return the second derivative of total RF voltage.

deltaVRF(z, I0)

Return the generator voltage minus beam loading voltage.

update_feedback()

Force feedback update from current CavityResonator parameters.

sample_voltage()

Sample the voltage seen by a zero charge particle during an RF period.

References

[1] Wilson, P. B. (1994). Fundamental-mode rf design in e+ e− storage ring factories. In Frontiers of Particle Beams: Factories with e+ e-Rings (pp. 293-311). Springer, Berlin, Heidelberg.

[2] Naoto Yamamoto, Alexis Gamelin, and Ryutaro Nagaoka. “Investigation of Longitudinal Beam Dynamics With Harmonic Cavities by Using the Code Mbtrack.” IPAC’19, Melbourne, Australia, 2019.

__init__(ring: Synchrotron, m: int | float, Rs: float, Q: float, QL: float, detune: float, Ncav: int = 1, Vc: float = 0, theta: float = 0, n_bin: int = 75)[source]
init_tracking(beam: Beam)[source]

Initialization of the tracking.

Parameters

beam : Beam object

track(beam: Beam)[source]

Track a Beam object through the CavityResonator object.

Can be used with or without mpi. If used with mpi, beam.mpi.share_distributions must be called before.

The beam phasor is given at t=0 (synchronous particle) of the first non empty bunch.

Parameters

beam : Beam object

init_phasor_track(beam: Beam)[source]

Initialize the beam phasor for a given beam distribution using a tracking like method.

Follow the same steps as the track method but in the “rf” reference frame and without any modifications on the beam.

Parameters

beam : Beam object

phasor_decay(time: float, ref_frame: str = 'beam')[source]

Compute the beam phasor decay during a given time span, assuming that no particles are crossing the cavity during the time span.

Parameters

timefloat

Time span in [s], can be positive or negative.

ref_framestring, optional

Reference frame to be used, can be “beam” or “rf”.

phasor_evol(profile: ndarray[tuple[int, ...], dtype[_ScalarType_co]], bin_length: float, charge_per_mp: float, ref_frame: str = 'beam')[source]

Compute the beam phasor evolution during the crossing of a bunch using an analytic formula [1].

Assume that the phasor decay happens before the beam loading.

Parameters

profilearray

Longitudinal profile of the bunch in [number of macro-particle].

bin_lengthfloat

Length of a bin in [s].

charge_per_mpfloat

Charge per macro-particle in [C].

ref_framestring, optional

Reference frame to be used, can be “beam” or “rf”.

References

[1] mbtrack2 manual.

init_phasor(beam: Beam)[source]

Initialize the beam phasor for a given beam distribution using an analytic formula [1].

No modifications on the Beam object.

Parameters

beam : Beam object

References

[1] mbtrack2 manual.

property generator_phasor: complex

Generator phasor in [V]

property cavity_phasor: complex

Cavity total phasor in [V]

property cavity_phasor_record: ndarray[tuple[int, ...], dtype[_ScalarType_co]]

Last cavity phasor value of each bunch in [V]

property ig_phasor_record: ndarray[tuple[int, ...], dtype[_ScalarType_co]]

Last current generator phasor of each bunch in [A]

property DFB_ig_phasor: ndarray[tuple[int, ...], dtype[_ScalarType_co]]

Last direct feedback current generator phasor of each bunch in [A]

property cavity_voltage: float

Cavity total voltage in [V]

property cavity_phase: float

Cavity total phase in [rad]

property beam_voltage: float

Beam loading voltage in [V]

property beam_phase: float

Beam loading phase in [rad]

property loss_factor: float

Cavity loss factor in [V/C]

property m: int | float

Harmonic number of the cavity

property Ncav: int

Number of cavities

property Rs_per_cavity: float

Shunt impedance of a single cavity in [Ohm], defined as 0.5*Vc*Vc/Pc.

property Rs: float

Shunt impedance [ohm]

property RL: float

Loaded shunt impedance [ohm]

property Q: float

Quality factor

property QL: float

Loaded quality factor

property beta: float

Coupling coefficient

property detune: float

Cavity detuning [Hz] - defined as (fr - m*f1)

property fr: float

Resonance frequency of the cavity in [Hz]

property wr: float

Angular resonance frequency in [Hz.rad]

property psi: float

Tuning angle in [rad]

property filling_time: float

Cavity filling time in [s]

property Pc: float

Power dissipated in the cavity walls in [W]

Pb(I0: float) float[source]

Return power transmitted to the beam in [W] - near Eq. (4.2.3) in [1].

Parameters

I0float

Beam current in [A].

Returns

float

Power transmitted to the beam in [W].

Pr(I0: float) float[source]

Power reflected back to the generator in [W].

Parameters

I0float

Beam current in [A].

Returns

float

Power reflected back to the generator in [W].

Vbr(I0: float) float[source]

Return beam voltage at resonance in [V].

Parameters

I0float

Beam current in [A].

Returns

float

Beam voltage at resonance in [V].

Vb(I0: float) float[source]

Return beam voltage in [V].

Parameters

I0float

Beam current in [A].

Returns

float

Beam voltage in [V].

Z(f: float) complex[source]

Cavity impedance in [Ohm] for a given frequency f in [Hz]

set_optimal_detune(I0: float, F: float = 1)[source]

Set detuning to optimal conditions - second Eq. (4.2.1) in [1].

Parameters

I0float

Beam current in [A].

set_optimal_coupling(I0: float)[source]

Set coupling to optimal value - Eq. (4.2.3) in [1].

Parameters

I0float

Beam current in [A].

set_generator(I0: float)[source]

Set generator parameters (Pg, Vgr, theta_gr, Vg and theta_g) for a given current and set of parameters.

Parameters

I0float

Beam current in [A].

plot_phasor(I0: float, ax: Axes | None = None) Axes[source]

Plot phasor diagram showing the vector addition of generator and beam loading voltage.

Parameters

I0float

Beam current in [A].

axAxes, optional

Polar Axes where the plot is displayed. If None, a new polar figure is created.

Returns

ax : Axes

is_CBI_stable(I0: float, synchrotron_tune: float | None = None, modes: list[float] | float | None = None, bool_return: bool = False) tuple[float, int] | bool[source]

Check Coupled-Bunch-Instability stability. Effect of Direct RF feedback is not included.

This method is a wraper around lcbi_growth_rate to caluclate the CBI growth rate from the CavityResonator own impedance.

See lcbi_growth_rate for details.

Parameters

I0float

Beam current in [A].

synchrotron_tunefloat, optinal

Fractional number of longitudinal tune. If None, synchrotron_tune is computed using self.Vc as total voltage. Default is None.

modesfloat or list of float, optional

Coupled Bunch Instability mode numbers to consider. If not None, return growth_rates of the input coupled bunch modes. Default is None.

bool_returnbool, optional
If True:
  • return True if the most unstable mode is stronger than

longitudinal damping time. - return False if it is not the case.

Default is False.

Returns

growth_ratefloat

Maximum coupled bunch instability growth rate in [s-1].

muint

Coupled bunch mode number corresponding to the maximum coupled bunch instability growth rate.

is_DC_Robinson_stable(I0: float) bool[source]

Check DC Robinson stability - Eq. (6.1.1) [1]

Parameters

I0float

Beam current in [A].

Returns

bool

plot_DC_Robinson_stability(detune_range: _Buffer | _SupportsArray[dtype[Any]] | _NestedSequence[_SupportsArray[dtype[Any]]] | bool | int | float | complex | str | bytes | _NestedSequence[bool | int | float | complex | str | bytes] = [-100000.0, 100000.0], ax: Axes | None = None) Axes[source]

Plot DC Robinson stability limit.

Parameters

detune_rangelist or array, optional

Range of tuning to plot in [Hz].

axAxes, optional

Axes where the plot is displayed. If None, a new figure is created.

Returns

ax : Axes

VRF(z: ndarray[tuple[int, ...], dtype[_ScalarType_co]], I0: float, F: float = 1, PHI: float = 0)[source]

Total RF voltage taking into account form factor amplitude F and form factor phase PHI

dVRF(z: ndarray[tuple[int, ...], dtype[_ScalarType_co]], I0: float, F: float = 1, PHI: float = 0)[source]

Return derivative of total RF voltage taking into account form factor amplitude F and form factor phase PHI

ddVRF(z: ndarray[tuple[int, ...], dtype[_ScalarType_co]], I0: float, F: float = 1, PHI: float = 0)[source]

Return the second derivative of total RF voltage taking into account form factor amplitude F and form factor phase PHI

deltaVRF(z: ndarray[tuple[int, ...], dtype[_ScalarType_co]], I0: float, F: float = 1, PHI: float = 0)[source]

Return the generator voltage minus beam loading voltage taking into account form factor amplitude F and form factor phase PHI

update_feedback()[source]

Force feedback update from current CavityResonator parameters.

sample_voltage(n_points: int | float = 10000.0, index: int = 0)[source]

Sample the voltage seen by a zero charge particle during an RF period.

Parameters

n_pointsint or float, optional

Number of sample points. The default is 1e4.

indexint, optional

RF bucket number to sample. Be carful if index > 0 as no new beam loading is not taken into account here. The default is 0.

Returns

posarray of float

Array of position from -T1/2 to T1/2.

voltage_recarray of float

Recoring of the voltage.

class ProportionalLoop(ring: Synchrotron, cav_res: CavityResonator, gain_A: float, gain_P: float, delay: int)[source]

Bases: object

Proportional feedback loop to control a CavityResonator amplitude and phase.

Feedback setpoints are cav_res.Vc and cav_res.theta.

The loop must be added to the CavityResonator object using:

cav_res.feedback.append(loop)

Parameters

ring : Synchrotron object cav_res : CavityResonator

The CavityResonator which is to be controlled.

gain_Afloat

Amplitude (voltage) gain of the feedback.

gain_Pfloat

Phase gain of the feedback.

delayint

Feedback delay in unit of turns. Must be supperior or equal to 1.

__init__(ring: Synchrotron, cav_res: CavityResonator, gain_A: float, gain_P: float, delay: int)[source]
track()[source]

Tracking method for the amplitude and phase loop.

Returns

None.

class TunerLoop(ring: Synchrotron, cav_res: CavityResonator, gain: float = 0.01, avering_period: int | None = None, offset: float = 0)[source]

Bases: object

Cavity tuner loop used to control a CavityResonator tuning (psi or detune) in order to keep the phase between cavity and generator current constant.

Only a proportional controller is assumed.

The loop must be added to the CavityResonator object using:

cav_res.feedback.append(loop)

Parameters

ring : Synchrotron object cav_res : CavityResonator

The CavityResonator which is to be controlled.

gainfloat

Proportional gain of the tuner loop. If not specified, 0.01 is used.

avering_periodint, optional

Period during which the phase difference is monitored and averaged. Then the feedback correction is applied every avering_period turn. Unit is turn number. A value longer than one synchrotron period (1/fs) is recommended. If None, 2-synchrotron period (2/fs) is used, although it is too fast compared to the actual situation. Default is None.

offsetfloat, optional

Tuning offset in [rad]. Default is 0.

__init__(ring: Synchrotron, cav_res: CavityResonator, gain: float = 0.01, avering_period: int | None = None, offset: float = 0)[source]
track()[source]

Tracking method for the tuner loop.

Returns

None.

class ProportionalIntegralLoop(ring: Synchrotron, cav_res: CavityResonator, gain: list[float], sample_num: int, every: int, delay: int, IIR_cutoff: float = 0, FF: bool = True)[source]

Bases: object

Proportional Integral (PI) loop to control a CavityResonator amplitude and phase via generator current (ig) [1,2].

Feedback reference targets (setpoints) are cav_res.Vc and cav_res.theta.

The ProportionalIntegralLoop should be initialized only after generator parameters are set.

The loop must be added to the CavityResonator object using:

cav_res.feedback.append(loop)

When the reference target is changed, it might be needed to re-initialize the feedforward constant by calling init_FFconst().

Parameters

ring : Synchrotron object cav_res : CavityResonator object

CavityResonator in which the loop will be added.

gainfloat list

Proportional gain (Pgain) and integral gain (Igain) of the feedback system in the form of a list [Pgain, Igain]. In case of normal conducting cavity (QL~1e4), the Pgain of ~1.0 and Igain of ~1e4(5) are usually used. In case of super conducting cavity (QL > 1e6), the Pgain of ~100 can be used. In a “bad” parameter set, unstable oscillation of the cavity voltage can be caused. So, a parameter scan of the gain should be made. Igain * ring.T1 / dtau is Ki defined as a the coefficients for integral part in of [1], where dtau is a clock period of the PI controller.

sample_numint

Number of bunch over which the mean cavity voltage is computed. Units are in bucket numbers.

everyint

Sampling and clock period of the feedback controller Time interval between two cavity voltage monitoring and feedback. Units are in bucket numbers.

delayint

Loop delay of the PI feedback system. Units are in bucket numbers.

IIR_cutofffloat, optinal

Cutoff frequency of the IIR filter in [Hz]. If 0, cutoff frequency is infinity. Default is 0.

FFbool, optional

Boolean switch to use feedforward constant. True is recommended to prevent a cavity voltage drop in the beginning of the tracking. In case of small Pgain (QL ~ 1e4), cavity voltage drop may cause beam loss due to static Robinson. Default is True.

Methods

track()

Tracking method for the Cavity PI control feedback.

init_Ig2Vg_matrix()

Initialize matrix for Ig2Vg_matrix.

init_FFconst()

Initialize feedforward constant.

Ig2Vg_matrix()

Return Vg from Ig using matrix formalism.

Ig2Vg()

Go from Ig to Vg and apply values.

Vg2Ig(Vg)

Return Ig from Vg (assuming constant Vg).

IIR_init(cutoff)

Initialization for the IIR filter.

IIR(input)

Return IIR filter output.

IIRcutoff()

Return IIR cutoff frequency in [Hz].

Notes

Assumption : delay >> every ~ sample_num

Adjusting ig(Vg) parameter to keep the cavity voltage constant (to target values). The “track” method is

1) Monitoring the cavity voltage phasor Mean cavity voltage value between specified bunch number (sample_num) is monitored with specified interval period (every). 2) Changing the ig phasor The ig phasor is changed according the difference of the specified reference values with specified gain (gain). By using ig instead of Vg, the cavity response can be taken account. 3) ig changes are reflected to Vg after the specifed delay (delay) of the system

Vc–>(rot)–>IIR–>(-)–>(V->I,fac)–>PI–>Ig –> Vg

Ref

Examples

PF; E0=2.5GeV, C=187m, frf=500MHz

QL=11800, fs0=23kHz ==> gain=[0.5,1e4], sample_num=8, every=7(13ns), delay=500(1us)

is one reasonable parameter set.

The practical clock period is 13ns.

==> Igain_PF = Igain_mbtrack * Trf / 13ns = Igain_mbtrack * 0.153

References

[1] : https://en.wikipedia.org/wiki/PID_controller [2] : Yamamoto, N., Takahashi, T., & Sakanaka, S. (2018). Reduction and compensation of the transient beam loading effect in a double rf system of synchrotron light sources. PRAB, 21(1), 012001.

__init__(ring: Synchrotron, cav_res: CavityResonator, gain: list[float], sample_num: int, every: int, delay: int, IIR_cutoff: float = 0, FF: bool = True)[source]
track(apply_changes: bool = True)[source]

Tracking method for the Cavity PI control feedback.

Returns

None.

init_Ig2Vg_matrix()[source]

Initialize matrix for Ig2Vg_matrix.

Shoud be called before first use of Ig2Vg_matrix and after each cavity parameter change.

init_FFconst()[source]

Initialize feedforward constant.

Ig2Vg_matrix() ndarray[tuple[int, ...], dtype[_ScalarType_co]][source]

Return Vg from Ig using matrix formalism. Warning: self.init_Ig2Vg should be called after each CavityResonator parameter change.

Ig2Vg()[source]

Go from Ig to Vg.

Apply new values to cav_res.generator_phasor_record, cav_res.Vg and cav_res.theta_g from ig_phasor_record.

Vg2Ig(Vg: float) float[source]

Return Ig from Vg (assuming constant Vg).

Eq.25 of ref [2] assuming the dVg/dt = 0.

IIR_init(cutoff: float)[source]

Initialization for the IIR filter.

Parameters

cutofffloat

Cutoff frequency of the IIR filter in [Hz]. If 0, cutoff frequency is infinity.

IIR(input: ndarray[tuple[int, ...], dtype[_ScalarType_co]]) ndarray[tuple[int, ...], dtype[_ScalarType_co]][source]

Return IIR filter output.

property IIRcutoff: float

Return IIR cutoff frequency in [Hz].

class DirectFeedback(DFB_gain: float, DFB_phase_shift: float, DFB_sample_num: int | None = None, DFB_every: int | None = None, DFB_delay: int | None = None, **kwargs)[source]

Bases: ProportionalIntegralLoop

Direct RF FB (DFB) via generator current (ig) using PI controller [1,2].

The DirectFeedback inherits from ProportionalIntegralLoop so all mandatory parameters from ProportionalIntegralLoop should be passed as kwargs when initializing a DirectFeedback object.

To avoid cavity-beam unmatching (large synchrotron oscilation of beam), the CavityResonator generator parameters should be set before initialization.

The loop must be added to the CavityResonator object using:

cav_res.feedback.append(loop)

Parameters

DFB_gainfloat

Gain of the DFB.

DFB_phase_shiftfloat

Phase shift of the DFB.

DFB_sample_numint, optional

Sample number to monitor Vc for the DFB. The averaged Vc value in DFB_sample_num is monitored. Units are in bucket numbers. If None, value from the ProportionalIntegralLoop is used. Default is None.

DFB_everyint, optional

Interval to monitor and change Vc for the DFB. Units are in bucket numbers. If None, value from the ProportionalIntegralLoop is used. Default is None.

DFB_delayint, optional

Loop delay of the DFB. Units are in bucket numbers. If None, value from the ProportionalIntegralLoop is used. Default is None.

Attributes

phase_shiftfloat

Phase shift applied to DFB, defined as psi - DFB_phase_shift.

DFB_psi: float

Return detuning angle value with Direct RF feedback in [rad].

DFB_alphafloat

Return the amplitude ratio alpha of the DFB.

DFB_gammafloat

Return the gain gamma of the DFB.

DFB_Rsfloat

Return the shunt impedance of the DFB in [ohm].

Methods

DFB_parameter_set(DFB_gain, DFB_phase_shift)

Set DFB gain and phase shift parameters.

track()

Tracking method for the Direct RF feedback.

DFB_Vg()

Return the generator voltage main and DFB components in [V].

DFB_fs()

Return the modified synchrotron frequency in [Hz].

References

[1] : Akai, K. (2022). Stability analysis of rf accelerating mode with feedback loops under heavy beam loading in SuperKEKB. PRAB, 25(10), 102002. [2] : N. Yamamoto et al. (2023). Stability survey of a double RF system with RF feedback loops for bunch lengthening in a low-emittance synchrotron ring. In Proc. IPAC’23. doi:10.18429/JACoW-IPAC2023-WEPL161

__init__(DFB_gain: float, DFB_phase_shift: float, DFB_sample_num: int | None = None, DFB_every: int | None = None, DFB_delay: int | None = None, **kwargs)[source]
property DFB_phase_shift: float

Return DFB phase shift.

property phase_shift: float

Return phase shift applied to DFB. Defined as self.cav_res.psi - self.DFB_phase_shift.

property DFB_psi: float

Return detuning angle value with Direct RF feedback in [rad].

Fig.4 of ref [1].

property DFB_alpha: float

Return the amplitude ratio alpha of the DFB.

Near Eq. (13) of [1].

property DFB_gamma: float

Return the gain gamma of the DFB.

Near Eq. (13) of [1].

property DFB_Rs: float

Return the shunt impedance of the DFB in [ohm].

Eq. (15) of [1].

DFB_parameter_set(DFB_gain: float, DFB_phase_shift: float)[source]

Set DFB gain and phase shift parameters.

Parameters

DFB_gainfloat

Gain of the DFB.

DFB_phase_shiftfloat

Phase shift of the DFB.

track()[source]

Tracking method for the Direct RF feedback.

Returns

None.

DFB_Vg(vc: float = -1) tuple[float, float][source]

Return the generator voltage main and DFB components in [V].

DFB_fs(vg_main: float = -1, vg_drf: float = -1) float[source]

Return the modified synchrotron frequency in [Hz].