mbtrack2.utilities.beamloading module¶
Module where the BeamLoadingEquilibrium class is defined.
- class BeamLoadingEquilibrium(ring: Synchrotron, cavity_list: list[CavityResonator], I0: float, auto_set_MC_theta: bool = False, F: list[float] | None = None, PHI: list[float] | None = None, B1: float = -0.2, B2: float = 0.2, N: int = 1000)[source]¶
Bases:
objectClass used to compute beam equilibrium profile for a given storage ring and a list of RF cavities of any harmonic. The class assumes an uniform filling of the storage ring. Based on [1] which is an extension of [2].
Parameters¶
ring : Synchrotron object cavity_list : list of CavityResonator objects I0 : beam current in [A]. auto_set_MC_theta : if True, allow class to change cavity phase for
CavityResonator objetcs with m = 1 (i.e. main cavities)
F : list of form factor amplitude PHI : list of form factor phase B1 : lower intergration boundary B2 : upper intergration boundary N : int or float, optinal
Number of points used for longitudinal grid. Default is 1000.
References¶
[1] Gamelin, A., & Yamamoto, N. (2021). Equilibrium bunch density distribution with multiple active and passive RF cavities. IPAC’21 (MOPAB069). [2] Venturini, M. (2018). Passive higher-harmonic rf cavities with general settings and multibunch instabilities in electron storage rings. Physical Review Accelerators and Beams, 21(11), 114404.
- __init__(ring: Synchrotron, cavity_list: list[CavityResonator], I0: float, auto_set_MC_theta: bool = False, F: list[float] | None = None, PHI: list[float] | None = None, B1: float = -0.2, B2: float = 0.2, N: int = 1000)[source]¶
- energy_balance() float[source]¶
Return energy balance for the synchronous particle (z = 0 ,delta = 0).
- u(z: ndarray[tuple[int, ...], dtype[_ScalarType_co]]) ndarray[tuple[int, ...], dtype[_ScalarType_co]][source]¶
Scaled potential u(z)
- du_dz(z: ndarray[tuple[int, ...], dtype[_ScalarType_co]]) ndarray[tuple[int, ...], dtype[_ScalarType_co]][source]¶
Partial derivative of the scaled potential u(z) by z
- uexp(z: ndarray[tuple[int, ...], dtype[_ScalarType_co]]) ndarray[tuple[int, ...], dtype[_ScalarType_co]][source]¶
- integrate_func(f: ndarray[tuple[int, ...], dtype[_ScalarType_co]], g: ndarray[tuple[int, ...], dtype[_ScalarType_co]]) ndarray[tuple[int, ...], dtype[_ScalarType_co]][source]¶
Return Integral[f*g]/Integral[f] between B1 and B2
- to_solve(x: ndarray[tuple[int, ...], dtype[_ScalarType_co]] | list, CM: bool = True) ndarray[tuple[int, ...], dtype[_ScalarType_co]][source]¶
System of non-linear equation to solve to find the form factors F and PHI at equilibrum. The system is composed of Eq. (B6) and (B7) of [1] for each cavity. If auto_set_MC_theta == True, the system also find the main cavity phase to impose energy balance or cancel center of mass offset. If CM is True, the system imposes zero center of mass offset, if False, the system imposes energy balance.
- voltage(z: ndarray[tuple[int, ...], dtype[_ScalarType_co]]) ndarray[tuple[int, ...], dtype[_ScalarType_co]][source]¶
Return the RF system total voltage at position z
- dV(z: ndarray[tuple[int, ...], dtype[_ScalarType_co]]) ndarray[tuple[int, ...], dtype[_ScalarType_co]][source]¶
Return derivative of the RF system total voltage at position z
- ddV(z: ndarray[tuple[int, ...], dtype[_ScalarType_co]]) ndarray[tuple[int, ...], dtype[_ScalarType_co]][source]¶
Return the second derivative of the RF system total voltage at position z
- deltaVRF(z: ndarray[tuple[int, ...], dtype[_ScalarType_co]]) ndarray[tuple[int, ...], dtype[_ScalarType_co]][source]¶
Return the generator voltage minus beam loading voltage of the total RF system at position z
- beam_equilibrium(x0: ndarray[tuple[int, ...], dtype[_ScalarType_co]] | list | None = None, tol: float = 0.0001, method: str = 'hybr', options: dict | None = None, plot: bool = False, CM: bool = True, verbiose: bool = False)[source]¶
Solve system of non-linear equation to find the form factors F and PHI at equilibrum. Can be used to compute the equilibrium bunch profile.
Parameters¶
x0 : initial guess tol : tolerance for termination of the algorithm method : method used by scipy.optimize.root to solve the system options : options given to scipy.optimize.root plot : if True, plot the equilibrium bunch profile CM : if True, the system imposes zero center of mass offset, if False, the system imposes energy balance. verbiose : if True, print out informations about convergence.
Returns¶
sol : OptimizeResult object representing the solution
- PTBL_threshold(I0: float, m: int | None = None, MC_index: int = 0, HC_index: int = 1) tuple[float, float, float][source]¶
Compute the periodic transient beam loading (PTBL) instability threshold based on [1].
The beam_equilibrium method should be called before the PTBL_threshold one.
Eq. (22) and (23) have been modified for cosine voltage convention for the main cavity.
Parameters¶
- I0float
Beam total current in [A].
- mint, optional
Number of bunches in the beam. The default is None. If None, uniform filling is assumed.
- MC_indexint, optional
Index of the main cavity in cavity_list. The default is 0.
- HC_indexint, optional
Index of the harmonic cavity in cavity_list. The default is 1.
Returns¶
- etafloat
Amplification factor, if bigger than 1 then the beam is possibly in the PTBL regime, Eq. (22) of [1].
- RQthfloat
R/Q threshold, Eq. (24) of [1].
- ffloat
f function, Eq. (23) of [1].
Reference¶
[1] : He, T. (2022). Novel perturbation method for judging the stability of the equilibrium solution in the presence of passive harmonic cavities. Physical Review Accelerators and Beams, 25(9), 094402.
- property R_factor: float¶
Touschek lifetime ratio as defined in [1].
Reference¶
[1] : Byrd, J. M., and M. Georgsson. “Lifetime increase using passive harmonic cavities in synchrotron light sources.” Physical Review Special Topics-Accelerators and Beams 4.3 (2001): 030701.
- property bunch_spectrum: ndarray[tuple[int, ...], dtype[_ScalarType_co]]¶
Return the bunch equilibrium spectrum.