mbtrack2.impedance.wakefield module

This module defines general classes to describes wakefields, impedances and wake functions.

class ComplexData(variable: ndarray[tuple[int, ...], dtype[_ScalarType_co]] = array([-1.e+15, 1.e+15]), function: ndarray[tuple[int, ...], dtype[_ScalarType_co]] = array([0, 0]))[source]

Bases: object

Define a general data structure for a complex function based on a pandas DataFrame.

Parameters

variablelist, numpy array

contains the function variable values

functionlist or numpy array of complex numbers

contains the values taken by the complex function

__init__(variable: ndarray[tuple[int, ...], dtype[_ScalarType_co]] = array([-1.e+15, 1.e+15]), function: ndarray[tuple[int, ...], dtype[_ScalarType_co]] = array([0, 0]))[source]
add(structure_to_add: ComplexData, method: str = 'zero', interp_kind: str = 'cubic', index_name: str = 'variable', sub: bool = False, rsub: bool = False) ComplexData[source]

Method to add two structures. If the data don’t have the same length, different cases are handled.

Parameters

structure_to_addComplexData object, int, float or complex.

Structure to be added.

methodstr, optional

Manage how the addition is done, possibilties are: -common: the common range of the index (variable) from the two

structures is used. The input data are cross-interpolated so that the result data contains the points at all initial points in the common domain.

-extrapolate: extrapolate the value of both ComplexData. -zero: outside of the common range, the missing values are zero. In

the common range, behave as “common” method.

The default is “zero”.

interp_kindstr, optional

Interpolation method which is passed to pandas and to scipy.interpolate.interp1d. The default is “cubic”.

index_namestr, optional

Name of the dataframe index passed to the method The default is “variable”.

subbool, optional

If True, substract structure_to_add from self. The default is False.

rsubbool, optional

If True, substract self from structure_to_add. The default is False.

Returns

ComplexData

Contains the sum of the two inputs.

multiply(factor: float | int) ComplexData[source]

Multiply a data strucure with a float or an int. If the multiplication is done with something else, throw a warning.

initialize_coefficient()[source]

Define the impedance coefficients and the plane of the impedance that are presents in attributes of the class.

static name_and_coefficients_table() DataFrame[source]

Return a table associating the human readbale names of an impedance component and its associated coefficients and plane.

property power_x: float
property power_y: float
property component_type: str
class WakeFunction(variable: ndarray[tuple[int, ...], dtype[_ScalarType_co]] = array([-1.e+15, 1.e+15]), function: ndarray[tuple[int, ...], dtype[_ScalarType_co]] = array([0, 0]), component_type: str = 'long')[source]

Bases: ComplexData

Define a WakeFunction object based on a ComplexData object.

Parameters

variablearray-like

Time domain structure of the wake function in [s].

functionarray-like

Wake function values in [V/C] for “long” and in [V/C/m] for others.

component_typestr, optinal

Type of the wake function: “long”, “xdip”, “xquad”, …

Attributes

data : DataFrame wake_type : str

Methods

to_impedance()

Return an Impedance object from the WakeFunction data.

deconvolution(freq_lim, sigma, mu)

Compute a wake function from wake potential data.

plot()

Plot the wake function data.

loss_factor(sigma)

Compute the loss factor or the kick factor for a Gaussian bunch.

__init__(variable: ndarray[tuple[int, ...], dtype[_ScalarType_co]] = array([-1.e+15, 1.e+15]), function: ndarray[tuple[int, ...], dtype[_ScalarType_co]] = array([0, 0]), component_type: str = 'long')[source]
add(structure_to_add: WakeFunction, method: str = 'zero', sub: bool = False, rsub: bool = False) WakeFunction[source]

Method to add two WakeFunction objects. The two structures are compared so that the addition is self-consistent.

multiply(factor: float | int) WakeFunction[source]

Multiply a WakeFunction object with a float or an int. If the multiplication is done with something else, throw a warning.

to_impedance(freq_lim: float, nout: int | float | None = None, sigma: float | None = None, mu: float | None = None) Impedance[source]

Return an Impedance object from the WakeFunction data. The WakeFunction data is assumed to be sampled equally.

If both sigma and mu are given, the impedance data is divided by the Fourier transform of a Gaussian distribution. It then can be used to go from wake potential data to impedance.

Parameters

freq_limfloat

The maximum frequency for calculating the impedance [Hz].

noutint or float, optional

Length of the transformed axis of the output. If None, it is taken to be 2*(n-1) where n is the length of the input. If nout > n, the input is padded with zeros. If nout < n, the input it cropped. Note that, for any nout value, nout//2+1 input points are required. The default is None.

sigmafloat, optional

RMS of the Gaussian distribution in [s]. The default is None.

mufloat, optional

Offset of the Gaussian distribution in [s]. The default is None.

deconvolution(freq_lim: float, sigma: float, mu: float, nout: int | float | None = None) WakeFunction[source]

Compute a wake function from wake potential data.

The present data loaded into the WakeFunction object is assumed to be an uniformly sampled wake potential from a Gaussian bunch distribution with RMS bunch length sigma and offset mu.

The offset mu depends on the code used to compute the wake potential:
  • for CST, it is the first time sample, so self.data.index[0].

  • for ABCI, it is “-1 * ISIG * SIG / c”.

Parameters

freq_limfloat

The maximum frequency for calculating the impedance [Hz].

sigmafloat

RMS of the Gaussian distribution in [s].

mufloat

Offset of the Gaussian distribution in [s].

noutint or float, optional

Length of the transformed axis of the output. If None, it is taken to be 2*(n-1) where n is the length of the input. If nout > n, the input is padded with zeros. If nout < n, the input it cropped. Note that, for any nout value, nout//2+1 input points are required. The default is None.

plot() Axes[source]

Plot the wake function data.

Returns

axclass:matplotlib.axes.Axes

Wake function data.

loss_factor(sigma: float) float[source]

Compute the loss factor or the kick factor for a Gaussian bunch. Formulas from Eq. (5.6), Eq. (5.12) and Eq. (5.17) of [1].

Parameters

sigmafloat

RMS bunch length in [s]

Returns

kloss: float

Loss factor in [V/C] or kick factor in [V/C/m] depanding on the impedance type.

References

[1] : Zotter, Bruno W., and Semyon A. Kheifets. Impedances and wakes in high-energy particle accelerators. World Scientific, 1998.

class Impedance(variable: ndarray[tuple[int, ...], dtype[_ScalarType_co]] = array([-1.e+15, 1.e+15]), function: ndarray[tuple[int, ...], dtype[_ScalarType_co]] = array([0, 0]), component_type: str = 'long')[source]

Bases: ComplexData

Define an Impedance object based on a ComplexData object.

Parameters

variablearray-like

Frequency domain structure of the impedance in [Hz].

functionarray-like

Impedance values in [Ohm].

component_typestr, optinal

Type of the impedance: “long”, “xdip”, “xquad”, …

Attributes

data : DataFrame impedance_type : str

Methods

loss_factor(sigma)

Compute the loss factor or the kick factor for a Gaussian bunch.

to_wakefunction()

Return a WakeFunction object from the impedance data.

plot()

Plot the impedance data.

__init__(variable: ndarray[tuple[int, ...], dtype[_ScalarType_co]] = array([-1.e+15, 1.e+15]), function: ndarray[tuple[int, ...], dtype[_ScalarType_co]] = array([0, 0]), component_type: str = 'long')[source]
add(structure_to_add: Impedance, beta_x: float = 1, beta_y: float = 1, method: str = 'zero', sub: bool = False, rsub: bool = False) Impedance[source]

Method to add two Impedance objects. The two structures are compared so that the addition is self-consistent. Beta functions can be precised as well.

multiply(factor: float | int) Impedance[source]

Multiply a Impedance object with a float or an int. If the multiplication is done with something else, throw a warning.

loss_factor(sigma: float) float[source]

Compute the loss factor or the kick factor for a Gaussian bunch. Formulas from Eq. (2) p239 and Eq.(7) p240 of [1].

Parameters

sigmafloat

RMS bunch length in [s]

Returns

kloss: float

Loss factor in [V/C] or kick factor in [V/C/m] depanding on the impedance type.

References

[1] : Handbook of accelerator physics and engineering, 3rd printing.

to_wakefunction(nout: int | float | None = None, trim: bool = False) WakeFunction[source]

Return a WakeFunction object from the impedance data. The impedance data is assumed to be sampled equally.

Parameters

noutint or float, optional

Length of the transformed axis of the output. If None, it is taken to be 2*(n-1) where n is the length of the input. If nout > n, the input is padded with zeros. If nout < n, the input it cropped. Note that, for any nout value, nout//2+1 input points are required.

trimbool or float, optional

If True, the pseudo wake function is trimmed at time=0 and is forced to zero where time<0. If False, the original result is preserved. If a float is given, the pseudo wake function is trimmed from time <= trim to 0.

plot() Axes[source]

Plot the impedance data.

Returns

axclass:matplotlib.axes.Axes

Impedance data.

class WakeField(structure_list: list[Impedance | WakeFunction] | None = None, name: str | None = None)[source]

Bases: object

Defines a WakeField which corresponds to a single physical element which produces different types of wakes, represented by their Impedance or WakeFunction objects.

Parameters

structure_listlist, optional

list of Impedance/WakeFunction objects to add to the WakeField.

name : str, optional

Attributes

impedance_componentsarray of str

Impedance components present for this element.

wake_componentsarray of str

WakeFunction components present for this element.

componentsarray of str

Impedance or WakeFunction components present for this element.

Methods

append_to_model(structure_to_add)

Add Impedance/WakeFunction to WakeField.

list_to_attr(structure_list)

Add list of Impedance/WakeFunction to WakeField.

add_wakefields(wake1, beta1, wake2, beta2)

Add two WakeFields taking into account beta functions.

add_several_wakefields(wakefields, beta)

Add a list of WakeFields taking into account beta functions.

drop(component)

Delete a component or a list of component from the WakeField.

save(file)

Save WakeField element to file.

load(file)

Load WakeField element from file.

save_to_pyheadtail(filename)

Save WakeField to PyHEADTAIL format.

__init__(structure_list: list[Impedance | WakeFunction] | None = None, name: str | None = None)[source]
_add(structure_to_add: WakeField) WakeField[source]

Allow to add two WakeField element with different components.

append_to_model(structure_to_add: Impedance | WakeFunction)[source]

Add Impedance/WakeFunction component to WakeField.

Parameters

structure_to_add : Impedance or WakeFunction object

list_to_attr(structure_list: list[Impedance | WakeFunction])[source]

Add list of Impedance/WakeFunction components to WakeField.

Parameters

structure_list : list of Impedance or WakeFunction objects.

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

Return an array of the impedance component names for the element.

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

Return an array of the wake function component names for the element.

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

Return an array of the component names for the element.

drop(component: str | list[str])[source]

Delete a component or a list of component from the WakeField.

Parameters

componentstr or list of str

Component or list of components to drop. If “Z”, drop all impedance components. If “W””, drop all wake function components.

Returns

None.

save(file: str)[source]

Save WakeField element to file.

The same pandas version is needed on both saving and loading computer for the pickle to work.

Parameters

filestr

File where the WakeField element is saved.

Returns

None.

save_to_pyheadtail(file: str)[source]

Saves wakefield model to pyheadtail format.

Parameters

filestr

filename to save the wakefield model

Returns

None.

static load(file: str) WakeField[source]

Load WakeField element from file.

The same pandas version is needed on both saving and loading computer for the pickle to work.

Parameters

filestr

File where the WakeField element is saved.

Returns

wakefieldWakeField

Loaded wakefield.

static add_wakefields(wake1: WakeField, beta1: ndarray[tuple[int, ...], dtype[_ScalarType_co]], wake2: WakeField, beta2: ndarray[tuple[int, ...], dtype[_ScalarType_co]]) WakeField[source]

Add two WakeFields taking into account beta functions.

Parameters

wake1WakeField

WakeField to add.

beta1array of shape (2,)

Beta functions at wake1 postion.

wake2WakeField

WakeField to add.

beta2array of shape (2,)

Beta functions at wake2 postion.

Returns

wake_sumWakeField

WakeField with summed components.

static add_several_wakefields(wakefields: list[WakeField], beta: ndarray[tuple[int, ...], dtype[_ScalarType_co]])[source]

Add a list of WakeFields taking into account beta functions.

Parameters

wakefieldslist of WakeField

WakeFields to add.

betaarray of shape (len(wakefields), 2)

Beta functions at the WakeField postions.

Returns

wake_sumWakeField

WakeField with summed components.