beam Package

beam Module

Module containing the fundamental beam class with methods to compute beam statistics

Authors

Danilo Quartullo, Helga Timko, ALexandre Lasheen

class blond.beam.beam.Beam(Ring, n_macroparticles, intensity)

Bases: object

Class containing the beam properties.

This class containes the beam coordinates (dt, dE) and the beam properties.

The beam coordinate ‘dt’ is defined as the particle arrival time to the RF station w.r.t. the reference time that is the sum of turns. The beam coordiate ‘dE’ is defined as the particle energy offset w.r.t. the energy of the synchronous particle.

The class creates a beam with zero dt and dE, see distributions to match a beam with respect to the RF and intensity effects.

Parameters
  • Ring (Ring) – Used to import different quantities such as the mass and the energy.

  • n_macroparticles (int) – total number of macroparticles.

  • intensity (float) – total intensity of the beam (in number of charge).

mass

mass of the particle [eV].

Type

float

charge

integer charge of the particle [e].

Type

int

beta

relativistic velocity factor [].

Type

float

gamma

relativistic mass factor [].

Type

float

energy

energy of the synchronous particle [eV].

Type

float

momentum

momentum of the synchronous particle [eV].

Type

float

dt

beam arrival times with respect to synchronous time [s].

Type

numpy_array, float

dE

beam energy offset with respect to the synchronous particle [eV].

Type

numpy_array, float

mean_dt

average beam arrival time [s].

Type

float

mean_dE

average beam energy offset [eV].

Type

float

sigma_dt

standard deviation of beam arrival time [s].

Type

float

sigma_dE

standard deviation of beam energy offset [eV].

Type

float

intensity

total intensity of the beam in number of charges [].

Type

float

n_macroparticles

total number of macroparticles in the beam [].

Type

int

ratio

ratio intensity per macroparticle [].

Type

float

n_macroparticles_lost

number of macro-particles marked as ‘lost’ [].

Type

int

id

unique macro-particle ID number; zero if particle is ‘lost’.

Type

numpy_array, int

See also

distributions.matched_from_line_density

match a beam with a given bunch profile.

distributions.matched_from_distribution_function

match a beam with a given distribution function in phase space.

Examples

>>> from input_parameters.ring import Ring
>>> from beam.beam import Beam
>>>
>>> n_turns = 10
>>> C = 100
>>> eta = 0.03
>>> momentum = 26e9
>>> ring = Ring(n_turns, C, eta, momentum, 'proton')
>>> n_macroparticle = 1e6
>>> intensity = 1e11
>>>
>>> my_beam = Beam(ring, n_macroparticle, intensity)
add_beam(other_beam)

Method to add the particles from another beam to this beam New particles are given id numbers sequential from last id of this beam Particles with id == 0 keep id == 0 and are included in addition

Parameters

other_beam (blond beam object) –

add_particles(new_particles)

Method to add array of new particles to beam object New particles are given id numbers sequential from last id of this beam

Parameters

new_particles (array-like) – (2, n) array of (dt, dE) for new particles

eliminate_lost_particles()

Eliminate lost particles from the beam coordinate arrays

gather(all=False)

MPI ONLY ROUTINE: Gather the beam coordinates to the master or all workers.

Parameters

all (boolean) – If true, every worker will get a copy of the whole beam coordinates. If false, only the master will gather the coordinates.

gather_losses(all=False)

MPI ONLY ROUTINE: Gather beam losses.

Parameters

all (boolean) – if true, all workers will gather the beam stats. If false, only the master will get the beam stats.

gather_statistics(all=False)

MPI ONLY ROUTINE: Gather beam statistics.

Parameters

all (boolean) – if true, all workers will gather the beam stats. If false, only the master will get the beam stats.

losses_below_energy(dE_min)

Beam losses based on lower energy cut.

Set to 0 all the particle’s id with dE below dE_min.

Parameters

dE_min (float) – minimum dE.

losses_energy_cut(dE_min, dE_max)

Beam losses based on energy cuts, e.g. on collimators.

Set to 0 all the particle’s id with dE not in the interval (dE_min, dE_max).

Parameters
  • dE_min (float) – minimum dE.

  • dE_max (float) – maximum dE.

losses_longitudinal_cut(dt_min, dt_max)

Beam losses based on longitudinal cuts.

Set to 0 all the particle’s id with dt not in the interval (dt_min, dt_max).

Parameters
  • dt_min (float) – minimum dt.

  • dt_max (float) – maximum dt.

losses_separatrix(Ring, RFStation)

Beam losses based on separatrix.

Set to 0 all the particle’s id not in the separatrix anymore.

Parameters
  • Ring (Ring) – Used to call the function is_in_separatrix.

  • RFStation (RFStation) – Used to call the function is_in_separatrix.

property n_macroparticles_alive

Number of transmitted macro-particles, defined as @property.

Returns

n_macroparticles_alive – number of macroparticles not lost.

Return type

int

property n_macroparticles_lost

Number of lost macro-particles, defined as @property.

Returns

n_macroparticles_lost – number of macroparticles lost.

Return type

int

split(random=False, fast=False)

MPI ONLY ROUTINE: Splits the beam equally among the workers for MPI processing.

Parameters
  • random (boolean) – Shuffle the beam before splitting, to be used with the approximation methonds.

  • fast (boolean) – If true, it assumes that every worker has already a copy of the beam so only the particle ids are distributed. If false, all the coordinates are distributed by the master to all the workers.

statistics()

Calculation of the mean and standard deviation of beam coordinates, as well as beam emittance using different definitions. Take no arguments, statistics stored in

  • mean_dt

  • mean_dE

  • sigma_dt

  • sigma_dE

class blond.beam.beam.Electron

Bases: blond.beam.beam.Particle

class blond.beam.beam.Particle(user_mass, user_charge)

Bases: object

class blond.beam.beam.Positron

Bases: blond.beam.beam.Particle

class blond.beam.beam.Proton

Bases: blond.beam.beam.Particle

coasting_beam Module

Module to generate coasting beam

Authors

Simon Albright

blond.beam.coasting_beam.generate_coasting_beam(Beam, t_start, t_stop, spread=0.001, spread_type='dp/p', energy_offset=0, distribution='gaussian', user_distribution=None, user_probability=None)

distributions Module

Module to generate distributions

Authors

Danilo Quartullo, Helga Timko, Alexandre Lasheen, Juan F. Esteban Mueller, Theodoros Argyropoulos, Joel Repond

blond.beam.distributions.X0_from_bunch_length(bunch_length, bunch_length_fit, X_grid, sorted_X_dE0, n_points_grid, time_potential_low_res, distribution_function_, distribution_type, distribution_exponent, beam, full_ring_and_RF)

Function to find the corresponding H0 or J0 for a given bunch length. Used by matched_from_distribution_function()

blond.beam.distributions.bigaussian(Ring, RFStation, Beam, sigma_dt, sigma_dE=None, seed=None, reinsertion=False)

Function generating a Gaussian beam both in time and energy coordinates. Fills Beam.dt and Beam.dE arrays.

Parameters
  • Ring (class) – A Ring type class

  • RFStation (class) – An RFStation type class

  • Beam (class) – A Beam type class

  • sigma_dt (float) – R.m.s. extension of the Gaussian in time

  • sigma_dE (float (optional)) – R.m.s. extension of the Gaussian in energy; default is None and will match the energy coordinate according to bucket height and sigma_dt

  • seed (int (optional)) – Fixed seed to have a reproducible distribution

  • reinsertion (bool (optional)) – Re-insert particles that are generated outside the separatrix into the bucket; default in False

blond.beam.distributions.distribution_function(action_array, dist_type, length, exponent=None)

Distribution function (formulas from Laclare).

blond.beam.distributions.line_density(coord_array, dist_type, bunch_length, bunch_position=0, exponent=None)

Line density

blond.beam.distributions.matched_from_distribution_function(beam, full_ring_and_RF, distribution_function_input=None, distribution_user_table=None, main_harmonic_option='lowest_freq', TotalInducedVoltage=None, n_iterations=1, n_points_potential=10000.0, n_points_grid=1000, dt_margin_percent=0.4, extraVoltageDict=None, seed=None, distribution_exponent=None, distribution_type=None, emittance=None, bunch_length=None, bunch_length_fit=None, distribution_variable='Hamiltonian', process_pot_well=True, turn_number=0)

Function to generate a beam by inputing the distribution function (by choosing the type of distribution and the emittance). The potential well is preprocessed to check for the min/max and center the frame around the separatrix. An error will be raised if there is not a full potential well (2 max and 1 min at least), or if there are several wells (more than 2 max and 1 min, this case will be treated in the future). An adjustable margin (40% by default) is applied in order to be able to catch the min/max of the potential well that might be on the edge of the frame. The slippage factor should be updated to take the higher orders. Outputs should be added in order for the user to check step by step if his bunch is going to be well generated. More detailed ‘step by step’ documentation should be implemented The user can input a custom distribution function by setting the parameter distribution_type = ‘user_input’ and passing the function in the parameter distribution_options[‘function’], with the following definition: distribution_function(action_array, dist_type, length, exponent=None). The user can also add an input table by setting the parameter distribution_type = ‘user_input_table’, distribution_options[‘user_table_action’] = array of action (in H or in J) and distribution_options[‘user_table_distribution’]

blond.beam.distributions.matched_from_line_density(beam, full_ring_and_RF, line_density_input=None, main_harmonic_option='lowest_freq', TotalInducedVoltage=None, plot=False, figdir='fig', half_option='first', extraVoltageDict=None, n_iterations=100, n_points_potential=10000.0, n_points_grid=1000, dt_margin_percent=0.4, n_points_abel=10000.0, bunch_length=None, line_density_type=None, line_density_exponent=None, seed=None, process_pot_well=True)

Function to generate a beam by inputing the line density. The distribution function is then reconstructed with the Abel transform and the particles randomly generated.

blond.beam.distributions.populate_bunch(beam, time_grid, deltaE_grid, density_grid, time_step, deltaE_step, seed)

Method to populate the bunch using a random number generator from the particle density in phase space.

distributions_multibunch Module

Module to generate multibunch distributions

Authors

Danilo Quartullo, Alexandre Lasheen, Theodoros Argyropoulos

blond.beam.distributions_multibunch.compute_H0(emittance, H, J)
blond.beam.distributions_multibunch.compute_X_grid(normalization_DeltaE, time_array, potential_well, distribution_variable)
blond.beam.distributions_multibunch.match_a_bunch(normalization_DeltaE, beam, potential_well_coordinates, potential_well, seed, distribution_options, full_ring_and_RF=None)
blond.beam.distributions_multibunch.match_beam_from_distribution(beam, FullRingAndRF, GeneralParameters, distribution_options, n_bunches, bunch_spacing_buckets, main_harmonic_option='lowest_freq', TotalInducedVoltage=None, n_iterations=1, n_points_potential=10000.0, dt_margin_percent=0.4, seed=None)

This function generates n equaly spaced bunches for a stationary distribution and try to match them with intensity effects.

The corresponding distributions are specified by their exponent:

\[g_0(J) \sim (1-J/J_0)^{\text{exponent}}}\]

Knowing the distribution, to generate the phase space: - Compute the potential U - The value of H can be computed thanks to U - The action J can be integrated over the whole phase space - 2piJ = emittance, this restrict the value of J0 (or H0) - with g0(H) we can randomize the macroparticles

blond.beam.distributions_multibunch.match_beam_from_distribution_multibatch(beam, FullRingAndRF, GeneralParameters, distribution_options, n_bunches, bunch_spacing_buckets, n_batch, batch_spacing_buckets, main_harmonic_option='lowest_freq', TotalInducedVoltage=None, n_iterations=1, n_points_potential=10000.0, dt_margin_percent=0.4, seed=None)

This function generates n equaly spaced bunches for a stationary distribution and try to match them with intensity effects.

Then it copies the batch n_batch times with spacing batch_spacing_buckets

The corresponding distributions are specified by their exponent:

\[g_0(J) \sim (1-J/J_0)^{\text{exponent}}}\]

Knowing the distribution, to generate the phase space: - Compute the potential U - The value of H can be computed thanks to U - The action J can be integrated over the whole phase space - 2piJ = emittance, this restrict the value of J0 (or H0) - with g0(H) we can randomize the macroparticles

blond.beam.distributions_multibunch.matched_from_distribution_density_multibunch(beam, Ring, FullRingAndRF, distribution_options_list, n_bunches, bunch_spacing_buckets, intensity_list=None, minimum_n_macroparticles=None, main_harmonic_option='lowest_freq', TotalInducedVoltage=None, n_iterations_input=1, plot_option=False, seed=None)

Function to generate a multi-bunch beam using the matched_from_distribution_density function for each bunch. The extra parameters to include are the number of bunches and the spacing between two bunches (assumed constant presently). Moreover, the distribution_options_list corresponds to the distribution_options of the matched_from_distribution_density function. It can be inputed as a dictionary just like the matched_from_distribution_density function (assuming the same parameters for all bunches), or as a list of length n_bunches to have different parameters for each bunch.

blond.beam.distributions_multibunch.matched_from_line_density_multibunch(beam, Ring, FullRingAndRF, line_density_options_list, n_bunches, bunch_spacing_buckets, intensity_list=None, minimum_n_macroparticles=None, main_harmonic_option='lowest_freq', TotalInducedVoltage=None, half_option='first', plot_option=False, seed=None)

Function to generate a multi-bunch beam using the matched_from_distribution_density function for each bunch. The extra parameters to include are the number of bunches and the spacing between two bunches (assumed constant presently). Moreover, the line_density_options_list corresponds to the distribution_options of the matched_from_line_density function. It can be inputed as a dictionary just like the matched_from_line_density function (assuming the same parameters for all bunches), or as a list of length n_bunches to have different parameters for each bunch.

profile Module

Module to compute the beam profile through slices

Authors

Danilo Quartullo, Alexandre Lasheen, Juan F. Esteban Mueller

class blond.beam.profile.CutOptions(cut_left=None, cut_right=None, n_slices=100, n_sigma=None, cuts_unit='s', RFSectionParameters=None)

Bases: object

This class groups all the parameters necessary to slice the phase space distribution according to the time axis, apart from the array collecting the profile which is defined in the constructor of the class Profile below.

Parameters
  • cut_left (float) – Left edge of the slicing (optional). A default value will be set if no value is given.

  • cut_right (float) – Right edge of the slicing (optional). A default value will be set if no value is given.

  • n_slices (int) – Optional input parameters, corresponding to the number of \(\sigma_{RMS}\) of the Beam to slice (this will overwrite any input of cut_left and cut_right).

  • n_sigma (float) – defines the left and right extremes of the profile in case those are not given explicitly

  • cuts_unit (str) – the unit of cut_left and cut_right, it can be seconds ‘s’ or radians ‘rad’

  • RFSectionParameters (object) – RFSectionParameters[0][0] is necessary for the conversion from radians to seconds if cuts_unit = ‘rad’. RFSectionParameters[0][0] is the value of omega_rf of the main harmonic at turn number 0

cut_left
Type

float

cut_right
Type

float

n_slices
Type

int

n_sigma
Type

float

cuts_unit
Type

str

RFSectionParameters
Type

object

edges

contains the edges of the slices

Type

float array

bin_centers

contains the centres of the slices

Type

float array

Examples

>>> from input_parameters.ring import Ring
>>> from input_parameters.rf_parameters import RFStation
>>> self.ring = Ring(n_turns = 1, ring_length = 100,
>>> alpha = 0.00001, momentum = 1e9)
>>> self.rf_params = RFStation(Ring=self.ring, n_rf=1, harmonic=[4620],
>>>                  voltage=[7e6], phi_rf_d=[0.])
>>> CutOptions = profileModule.CutOptions(cut_left=0, cut_right=2*np.pi,
>>> n_slices = 100, cuts_unit='rad', RFSectionParameters=self.rf_params)
convert_coordinates(value, input_unit_type)

Method to convert a value from ‘rad’ to ‘s’.

get_slices_parameters()

Reuturn all the computed parameters.

set_cuts(Beam=None)

Method to set self.cut_left, self.cut_right, self.edges and self.bin_centers attributes. The frame is defined by \(n\sigma_{RMS}\) or manually by the user. If not, a default frame consisting of taking the whole bunch +5% of the maximum distance between two particles in the bunch will be taken in each side of the frame.

track_cuts(Beam)

Track the slice frame (limits and slice position) as the mean of the bunch moves. Requires Beam statistics! Method to be refined!

class blond.beam.profile.FilterOptions(filterMethod=None, filterExtraOptions=None)

Bases: object

This class defines the filter to be used turn after turn to smooth the bunch profile.

Parameters
  • filterMethod (string) – The only option available is ‘chebishev’

  • filterExtraOptions (dictionary) – Parameters for the Chebishev filter (see the method beam_profile_filter_chebyshev in filters_and_fitting.py in the toolbox package)

filterMethod
Type

string

filterExtraOptions
Type

dictionary

class blond.beam.profile.FitOptions(fit_option=None, fitExtraOptions=None)

Bases: object

This class defines the method to be used turn after turn to obtain the position and length of the bunch profile.

Parameters
  • fit_method (string) – Current options are ‘gaussian’, ‘fwhm’ (full-width-half-maximum converted to 4 sigma gaussian bunch) and ‘rms’. The methods ‘gaussian’ and ‘rms’ give both 4 sigma.

  • fitExtraOptions (unknown) – For the moment no options can be passed into fitExtraOptions

fit_method
Type

string

fitExtraOptions
Type

unknown

class blond.beam.profile.OtherSlicesOptions(smooth=False, direct_slicing=False)

Bases: object

This class groups all the remaining options for the Profile class.

Parameters
  • smooth (boolean) – If set True, this method slices the bunch not in the standard way (fixed one slice all the macroparticles contribute with +1 or 0 depending if they are inside or not). The method assigns to each macroparticle a real value between 0 and +1 depending on its time coordinate. This method can be considered a filter able to smooth the profile.

  • direct_slicing (boolean) – If set True, the profile is calculated when the Profile class below is created. If False the user has to manually track the Profile object in the main file after its creation

smooth
Type

boolean

direct_slicing
Type

boolean

class blond.beam.profile.Profile(Beam, CutOptions=<blond.beam.profile.CutOptions object>, FitOptions=<blond.beam.profile.FitOptions object>, FilterOptions=<blond.beam.profile.FilterOptions object>, OtherSlicesOptions=<blond.beam.profile.OtherSlicesOptions object>)

Bases: object

Contains the beam profile and related quantities including beam spectrum, profile derivative.

Parameters
  • Beam (object) – Beam from which the profile has to be calculated

  • CutOptions (object) – Options for profile cutting (see above)

  • FitOptions (object) – Options to get profile position and length (see above)

  • FilterOptions (object) – Options to set a filter (see above)

  • OtherSlicesOptions (object) – All remaining options, like smooth histogram and direct slicing (see above)

Beam
Type

object

n_slices

number of slices to be used

Type

int

cut_left

left extreme of the profile

Type

float

cut_right

right extreme of the profile

Type

float

n_sigma

defines the left and right extremes of the profile in case those are not given explicitly

Type

float

edges

contains the edges of the slices

Type

float array

bin_centers

contains the centres of the slices

Type

float array

bin_size

lenght of one bin (or slice)

Type

float

n_macroparticles

contains the histogram (or profile); its elements are real if the smooth histogram tracking is used

Type

float array

beam_spectrum

contains the spectrum of the beam (arb. units)

Type

float array

beam_spectrum_freq

contains the frequencies on which the spectrum is computed [Hz]

Type

float array

operations

contains all the methods to be called every turn, like slice track, fitting, filtering etc.

Type

list

bunchPosition

profile position [s]

Type

float

bunchLength

profile length [s]

Type

float

filterExtraOptions
Type

unknown (see above)

Examples

>>> n_slices = 100
>>> CutOptions = profileModule.CutOptions(cut_left=0,
>>>       cut_right=self.ring.t_rev[0], n_slices = n_slices, cuts_unit='s')
>>> FitOptions = profileModule.FitOptions(fit_option='gaussian',
>>>                                        fitExtraOptions=None)
>>> filter_option = {'pass_frequency':1e7,
>>>    'stop_frequency':1e8, 'gain_pass':1, 'gain_stop':2,
>>>    'transfer_function_plot':False}
>>> FilterOptions = profileModule.FilterOptions(filterMethod='chebishev',
>>>         filterExtraOptions=filter_option)
>>> OtherSlicesOptions = profileModule.OtherSlicesOptions(smooth=False,
>>>                             direct_slicing = True)
>>> self.profile4 = profileModule.Profile(my_beam, CutOptions = CutOptions,
>>>                     FitOptions= FitOptions,
>>>                     FilterOptions=FilterOptions,
>>>                     OtherSlicesOptions = OtherSlicesOptions)
apply_filter()

It applies Chebishev filter to the profile.

apply_fit()

It applies Gaussian fit to the profile.

beam_profile_derivative(mode='gradient')

The input is one of the three available methods for differentiating a function. The two outputs are the bin centres and the discrete derivative of the Beam profile respectively.*

beam_spectrum_freq_generation(n_sampling_fft)

Frequency array of the beam spectrum

beam_spectrum_generation(n_sampling_fft)

Beam spectrum calculation

fwhm(shift=0)

Computation of the bunch length and position from the FWHM assuming Gaussian line density.

fwhm_multibunch(n_bunches, bunch_spacing_buckets, bucket_size_tau, bucket_tolerance=0.4, shift=0)

Computation of the bunch length and position from the FWHM assuming Gaussian line density for multibunch case.

reduce_histo(dtype=<class 'numpy.uint32'>)
rms()

Computation of the RMS bunch length and position from the line density (bunch length = 4sigma).

rms_multibunch(n_bunches, bunch_spacing_buckets, bucket_size_tau, bucket_tolerance=0.4)

Computation of the bunch length (4sigma) and position from RMS.

scale_histo()
set_slices_parameters()
track()

Track method in order to update the slicing along with the tracker. The kwargs are currently only needed to forward the reduce kw argument needed for the MPI version.

sparse_slices Module

Module to compute beam slicing for a sparse beam Only valid for cases with constant revolution and RF frequencies

Authors

Juan F. Esteban Mueller

class blond.beam.sparse_slices.SparseSlices(RFStation, Beam, n_slices, filling_pattern, tracker='C', direct_slicing=False)

Bases: object

This class instantiates a Slice object for each filled bucket according to the provided filling pattern. Each slice object will be of the size of an RF bucket and will have the same number of slices.

Beam

Import (reference) Beam

RFParams

Import (reference) RFStation

n_filled_buckets

Number of buckets to be sliced

n_slices

Number of slices per bucket

set_cuts()

Method to set the self.cut_left_array and self.cut_right_array properties, with the limits being an RF period. This is done as a pre-processing.