ocelot package

Subpackages

Module contents

general ocelot description

class ocelot.Aperture(xmax=inf, ymax=inf, dx=0, dy=0, type='rect', eid=None)

Bases: ocelot.cpbd.elements.Element

xmax - half size in horizontal plane in [m], ymax - half size in vertical plane in [m], type - “rect” or “elliptical”.

class ocelot.Beam(x=0, xp=0, y=0, yp=0)

Bases: object

property dg
property emit_xn
property emit_yn
property g
len()
property p
properties = ['g', 'dg', 'emit_xn', 'emit_yn', 'p', 'pz', 'px', 'py']
property px
property py
property pz
sizes()
to_array(nslice=100, window_len=None)
class ocelot.BeamTransform(tws=None, x_opt=None, y_opt=None)

Bases: ocelot.cpbd.physics_proc.PhysProc

Beam matching

apply(p_array, dz)

the method is called on every step.

Parameters
  • p_array

  • dz

Returns

property twiss
class ocelot.Bend(l=0.0, angle=0.0, k1=0.0, k2=0.0, e1=0.0, e2=0.0, tilt=0.0, gap=0.0, h_pole1=0.0, h_pole2=0.0, fint=0.0, fintx=None, eid=None)

Bases: ocelot.cpbd.elements.Element

bending magnet l - length of magnet in [m], angle - angle of bend in [rad], k1 - strength of quadrupole lens in [1/m^2], k2 - strength of sextupole lens in [1/m^3], tilt - tilt of lens in [rad], e1 - the angle of inclination of the entrance face [rad], e2 - the angle of inclination of the exit face [rad]. fint - fringe field integral fintx - allows (fintx > 0) to set fint at the element exit different from its entry value. gap - the magnet gap [m], NOTE in MAD and ELEGANT: HGAP = gap/2 h_pole1 - the curvature (1/r) of the entrance face h_pole1 - the curvature (1/r) of the exit face

class ocelot.CSR

Bases: ocelot.cpbd.physics_proc.PhysProc

coherent synchrotron radiation Attributes:

self.step = 1 [in Navigator.unit_step] - step of the CSR kick applying for beam (ParticleArray) self.sigma_min = 1.e-4 - minimal sigma if gauss filtering applied self.traj_step = 0.0002 [m] - trajectory step or, other words, integration step for calculation of the CSR-wake self.apply_step = 0.0005 [m] - step of the calculation CSR kick, to calculate average CSR kick

CSR_K1(i, traj, NdW, gamma=None)
Parameters
  • i – index of the trajectories points for the convolution kernel is calculated;

  • traj – trajectory. traj[0,:] - longitudinal coordinate, traj[1,:], traj[2,:], traj[3,:] - rectangular coordinates, traj[4,:], traj[5,:], traj[6,:] - tangential unit vectors

  • NdW – list [N, dW], NdW[0] is number of mesh points, NdW[1] = dW > 0 is increment. Mesh = Mesh = (-N:0) * dW

  • gamma

Returns

K0_fin_inf(i, traj, w_range, gamma)

Radiative interaction is coming from the infinite straight line that is assumed before the specified CSR region and it is calculated analytically

Parameters
  • i

  • traj

  • w_range

  • gamma

Returns

K0_inf_anf(i, traj, wmin)
Parameters
  • i – index of the trajectories points for the convolution kernel is calculated;

  • traj – trajectory. traj[0,:] - longitudinal coordinate, traj[1,:], traj[2,:], traj[3,:] - rectangular coordinates, traj[4,:], traj[5,:], traj[6,:] - tangential unit vectors

  • wmin – the first coordinate of the mash

Returns

K0_inf_inf(i, traj, w_range)

Radiative interaction is coming from the infinite straight line that is assumed before the specified CSR region and it is calculated analytically

Parameters
  • i

  • traj

  • w_range

Returns

apply(p_array, delta_s)

the method is called on every step.

Parameters
  • p_array

  • dz

Returns

finalize(*args, **kwargs)

the method is called at the end of tracking

Returns

plot_wake(p_array, lam_K1, itr_ra, s1, st)

Method to plot CSR wakes on each step and save pictures in the working folder. Might be time-consuming.

Parameters
  • p_array

  • lam_K1

  • itr_ra

  • s1

  • st

Returns

prepare(lat)

calculation of trajectory in rectangular coordinates calculation of the z_csr_start :param lat: Magnetic Lattice :return: self.csr_traj: trajectory. traj[0,:] - longitudinal coordinate,

traj[1,:], traj[2,:], traj[3,:] - rectangular coordinates, traj[4,:], traj[5,:], traj[6,:] - tangential unit vectors

class ocelot.Cavity(l=0.0, v=0.0, phi=0.0, freq=0.0, vx_up=0, vy_up=0, vxx_up=0, vxy_up=0, vx_down=0, vy_down=0, vxx_down=0, vxy_down=0, eid=None)

Bases: ocelot.cpbd.elements.Element

Standing wave RF cavity v - voltage [GV] freq - frequency [Hz] phi - phase in [deg] vx_{up/down}, vy_{up/down} - zero order kick of a {up/down}stream coupler vxx_{up/down}, vxy_{up/down} - first order kick a {up/down}stream coupler

class ocelot.CavityTM(v=0, freq=0.0, phi=0.0)

Bases: ocelot.cpbd.optics.TransferMap

map4cav(X, E, V, freq, phi, z=0)
class ocelot.Drift(l=0.0, eid=None)

Bases: ocelot.cpbd.elements.Element

drift - free space l - length of drift in [m]

class ocelot.EbeamParams(lattice, tws0, coupling=0.01, nsuperperiod=1)

Bases: object

integrals_id()
class ocelot.Edge(l=0.0, angle=0.0, k1=0.0, edge=0.0, tilt=0.0, dtilt=0.0, dx=0.0, dy=0.0, h_pole=0.0, gap=0.0, fint=0.0, pos=1, eid=None)

Bases: ocelot.cpbd.elements.Bend

class ocelot.Element(eid=None)

Bases: object

Element is a basic beamline building element Accelerator optics elements are subclasses of Element Arbitrary set of additional parameters can be attached if necessary

class ocelot.EmptyProc(step=1)

Bases: ocelot.cpbd.physics_proc.PhysProc

class ocelot.Hcor(l=0.0, angle=0.0, eid=None)

Bases: ocelot.cpbd.elements.RBend

horizontal corrector, l - length of magnet in [m], angle - angle of bend in [rad],

class ocelot.KickTM(angle=0.0, k1=0.0, k2=0.0, k3=0.0, nkick=1)

Bases: ocelot.cpbd.optics.TransferMap

kick(X, l, angle, k1, k2, k3, energy, nkick=1)

does not work for dipole

kick_apply(X, l, angle, k1, k2, k3, energy, nkick, dx, dy, tilt)
class ocelot.LSC(step=1)

Bases: ocelot.cpbd.physics_proc.PhysProc

Longitudinal Space Charge smooth_param - 0.1 smoothing parameter, resolution = np.std(p_array.tau())*smooth_param

apply(p_array, dz)

wakes in V/pC

Parameters
  • p_array

  • dz

Returns

imp_lsc(gamma, sigma, w, dz)

gamma - energy sigma - transverse RMS size of the beam w - omega = 2*pi*f

imp_step_lsc(gamma, rb, w, dz)

longitudinal space-charge impedance in case of a stepped profile bunch

gamma - energy rb - transverse radius of the beam w - omega = 2*pi*f

impedance2wake(f, y)

Fourier transform with exp(-iwt) f - Hz y - Om s - Meter w - V/C

wake2impedance(s, w)

Fourier transform with exp(iwt) s - Meter w - V/C f - Hz y - Om

wake_lsc(s, bunch, gamma, sigma, dz)
class ocelot.LaserHeater(step=1)

Bases: ocelot.cpbd.physics_proc.LaserModulator

class ocelot.LaserModulator(step=1)

Bases: ocelot.cpbd.physics_proc.PhysProc

apply(p_array, dz)

the method is called on every step.

Parameters
  • p_array

  • dz

Returns

lambda_ph(energy)

Wavelength of the laser pulse

Parameters

energy – in [GeV] - beam energy

Returns

wavelength in [m]

r56(energy)

Method calculate R56 of the undulator

Parameters

energy – in [GeV] - beam energy

Returns

R56 in [m]

class ocelot.MagneticLattice(sequence, start=None, stop=None, method=<ocelot.cpbd.optics.MethodTM object>)

Bases: object

sequence - list of the elements, start - first element of lattice. If None, then lattice starts from first element of sequence, stop - last element of lattice. If None, then lattice stops by last element of sequence, method = MethodTM() - method of the tracking.

find_indices(element)
get_sequence_part(start, stop)
update_endings(lat_index, element, body_elements, element_util)
update_transfer_maps()
class ocelot.Marker(eid=None)

Bases: ocelot.cpbd.elements.Element

class ocelot.Matrix(l=0.0, delta_e=0, eid=None, **kwargs)

Bases: ocelot.cpbd.elements.Element

Matrix element

l = 0 - m, length of the matrix element r = np.zeros((6, 6)) - R - elements, first order t = np.zeros((6, 6, 6)) - T - elements, second order delta_e = 0 - GeV, energy gain along the matrix element

class ocelot.MethodTM(params=None)

Bases: object

The class creates a transfer map for elements that depend on user-defined parameters (“parameters”). By default, the parameters = {“global”: TransferMap}, which means that all elements will have linear transfer maps. You can also specify different transfer maps for any type of element.

# use linear matrices for all elements except Sextupole which will have nonlinear kick map (KickTM) method = MethodTM() method.global_method = TransferMap method.params[Sextupole] = KickTM

# All elements are assigned matrices of the second order. # For elements for which there are no matrices of the second order are assigned default matrices, e.g. linear matrices. method2 = MethodTM() method2.global_method = SecondTM

create_tm(element)
set_tm(element, method)
class ocelot.Monitor(l=0.0, eid=None)

Bases: ocelot.cpbd.elements.Element

class ocelot.Multipole(kn=0.0, eid=None)

Bases: ocelot.cpbd.elements.Element

kn - list of strengths

class ocelot.Navigator(lattice)

Bases: object

Navigator defines step (dz) of tracking and which physical process will be applied during each step. lattice - MagneticLattice Attributes:

unit_step = 1 [m] - unit step for all physics processes

Methods:
add_physics_proc(physics_proc, elem1, elem2)

physics_proc - physics process, can be CSR, SpaceCharge or Wake, elem1 and elem2 - first and last elements between which the physics process will be applied.

activate_apertures(start=None, stop=None)

activate apertures if thea exist in the lattice from

Parameters
  • start – element, activate apertures starting form element ‘start’ element

  • stop – element, activate apertures up to ‘stop’ element

Returns

add_physics_proc(physics_proc, elem1, elem2)

Method adds Physics Process.

Parameters
  • physics_proc – PhysicsProc, e.g. SpaceCharge, CSR, Wake …

  • elem1 – the element in the lattice where to start applying the physical process.

  • elem2 – the element in the lattice where to stop applying the physical process, can be the same as starting element.

Returns

check_overjump(dz, processes, phys_steps)
check_proc_bounds(dz, proc_list, phys_steps, active_process)
get_next()
get_phys_procs()

method return list of all physics processes which were added

Returns

list, list of PhysProc(s)

get_proc_list()
go_to_start()
hard_edge_step(dz)
remove_used_processes(processes)

in case physics processes are applied and do not more needed they are removed from table

Parameters

processes – list of processes are about to apply

Returns

None

reset_position()

method to reset Navigator position. :return:

class ocelot.Octupole(l=0.0, k3=0.0, tilt=0.0, eid=None)

Bases: ocelot.cpbd.elements.Element

octupole k3 - strength of octupole lens in [1/m^4], l - length of lens in [m].

class ocelot.Particle(x=0.0, y=0.0, px=0.0, py=0.0, s=0.0, p=0.0, tau=0.0, E=0.0, q=0.0)

Bases: object

particle to be used for tracking

class ocelot.ParticleArray(n=0)

Bases: object

array of particles of fixed size; for optimized performance (x, x’ = px/p0),(y, y’ = py/p0),(ds = c*tau, p = dE/(p0*c)) p0 - momentum

class LostParticleRecorder(n)

Bases: object

Stores information about particles that are getting lost.

Attributes:

lost_particles: List of indices of deleted particle. Notes: The indices of the initial ParticleArray are used. lp_to_pos_hist: Histogram of number of lost particles to position s. [(pos, num of lost particles), …]

add(inds, position)
array2ex_list(p_list)
array2list()
delete_particles(inds, record=True)

Deletes particles from the particle array via index. :param inds: Indices that will be removed from the particle array. :param record: If record is true the deleted particles will be saved in self.lost_particle_recorder. :return:

list2array(p_list)
property n
p()
px()
py()
rescale2energy(energy)

Method to rescale beam coordinates with new energy

Parameters

energy – new energy

Returns

rm_particle(index)

Method removes a “bad” particle with particular “index”. :param index: :return:

rm_tails(xlim, ylim, px_lim, py_lim)

Method removes particles outside range [-xlim, +xlim], [-px_lim, +px_lim] …

size()
property t
tau()
thin_out(nth=10, n0=0)

Method to thin out the particle array in n-th times. Means every n-th particle will be saved in new Particle array

Parameters
  • nth – 10, every n-th particle will be taken to new Particle array

  • n0 – start from n0 particle

Returns

New ParticleArray

x()
y()
class ocelot.PhaseSpaceAperture(step=1)

Bases: ocelot.cpbd.physics_proc.PhysProc

Method to cut beam in longitudinal (by default), horizontal or/and vertical direction

Parameters
  • longitudinal – True, cutting in longitudinal direction

  • vertical – False, cutting in vertical direction

  • horizontal – False, cutting in horizontal direction

  • taumin – -5 longitudinal plane in [rms] from center of mass

  • taumax – 5 longitudinal plane in [rms] from center of mass

  • xmin – -5 horizontal plane in [rms] from center of mass

  • xmax – 5 horizontal plane in [rms] from center of mass

  • ymin – -5 vertical plane in [rms] from center of mass

  • ymax – 5 vertical plane in [rms] from center of mass

apply(p_array, dz)

the method is called on every step.

Parameters
  • p_array

  • dz

Returns

class ocelot.PhysProc(step=1)

Bases: object

Parent class for all Physics processes

Method prepare(self, lat)
  • the method is called at the moment of Physics Process addition to Navigator class.

Method apply(self, p_array, dz)
  • the method is called on every step.

Attribute step
  • number of steps in [Navigator.unit_step] self.step*Navigator.unit_step = [m]

Attribute indx0
  • number of start element in lattice.sequence - assigned in navigator.add_physics_proc()

Attribute indx1
  • number of stop element in lattice.sequence - assigned in navigator.add_physics_proc()

Attribute s_start
  • position of start element in lattice - assigned in navigator.add_physics_proc()

Attribute s_stop
  • position of stop element in lattice.sequence - assigned in navigator.add_physics_proc()

Attribute z0
  • current position of navigator - assigned in track.track() before p.apply()

apply(p_array, dz)

the method is called on every step.

Parameters
  • p_array

  • dz

Returns

finalize(*args, **kwargs)

the method is called at the end of tracking

Returns

prepare(lat)

method is called at the moment of Physics Process addition to Navigator class.

Parameters

lat

Returns

class ocelot.Quadrupole(l=0.0, k1=0, k2=0.0, tilt=0.0, eid=None)

Bases: ocelot.cpbd.elements.Element

quadrupole l - length of lens in [m], k1 - strength of quadrupole lens in [1/m^2], k2 - strength of sextupole lens in [1/m^3], tilt - tilt of lens in [rad].

class ocelot.RBend(l=0.0, angle=0.0, k1=0.0, k2=0.0, e1=None, e2=None, tilt=0.0, gap=0, h_pole1=0.0, h_pole2=0.0, fint=0.0, fintx=None, eid=None)

Bases: ocelot.cpbd.elements.Bend

rectangular bending magnet, l - length of magnet in [m], angle - angle of bend in [rad], k1 - strength of quadrupole lens in [1/m^2], k2 - strength of sextupole lens in [1/m^3], tilt - tilt of lens in [rad], e1 - the angle of inclination of the entrance face [rad], e2 - the angle of inclination of the exit face [rad]. fint - fringe field integral fintx - allows (fintx > 0) to set fint at the element exit different from its entry value. gap - the magnet gap [m], NOTE in MAD and ELEGANT: HGAP = gap/2 h_pole1 - the curvature (1/r) of the entrance face h_pole1 - the curvature (1/r) of the exit face

class ocelot.SBend(l=0.0, angle=0.0, k1=0.0, k2=0.0, e1=0.0, e2=0.0, tilt=0.0, gap=0, h_pole1=0.0, h_pole2=0.0, fint=0.0, fintx=None, eid=None)

Bases: ocelot.cpbd.elements.Bend

sector bending magnet, l - length of magnet in [m], angle - angle of bend in [rad], k1 - strength of quadrupole lens in [1/m^2], k2 - strength of sextupole lens in [1/m^3], tilt - tilt of lens in [rad], e1 - the angle of inclination of the entrance face [rad], e2 - the angle of inclination of the exit face [rad]. fint - fringe field integral fintx - allows (fintx > 0) to set fint at the element exit different from its entry value. gap - the magnet gap [m], NOTE in MAD and ELEGANT: HGAP = gap/2 h_pole1 - the curvature (1/r) of the entrance face h_pole1 - the curvature (1/r) of the exit face

class ocelot.SecondTM(r_z_no_tilt, t_mat_z_e)

Bases: ocelot.cpbd.optics.TransferMap

t_apply(R, T, X, dx, dy, tilt, U5666=0.0)
class ocelot.Sequence(l=0, refer=0)

Bases: object

class ocelot.Sextupole(l=0.0, k2=0.0, tilt=0.0, eid=None)

Bases: ocelot.cpbd.elements.Element

sextupole l - length of lens in [m], k2 - strength of sextupole lens in [1/m^3].

class ocelot.SmoothBeam

Bases: ocelot.cpbd.physics_proc.PhysProc

Physics Process for the beam smoothing. Can be applied when number of particles is not enough.

Attribute mslice

number of particles in the slice

# lat is the MagneticLattice navi = Navigator(lat)

smooth = SmoothBeam() smooth.mslice = 10000

navi.add_physics_process(smooth, start=elem, stop=elem) # elem is the lattice element where you want to apply smoothing

apply(p_array, dz)

the method is called on every step.

Parameters
  • p_array

  • dz

Returns

class ocelot.Solenoid(l=0.0, k=0.0, eid=None)

Bases: ocelot.cpbd.elements.Element

l - length in m, k - strength B0/(2B*rho)

class ocelot.SpaceCharge(step=1)

Bases: ocelot.cpbd.physics_proc.PhysProc

Space Charge physics process

Attributes:

self.step = 1 [in Navigator.unit_step] - step of the Space Charge kick applying self.nmesh_xyz = [63, 63, 63] - 3D mesh

Description:

The space charge forces are calculated by solving the Poisson equation in the bunch frame.

Then the Lorentz transformed electromagnetic field is applied as a kick in the laboratory frame. For the solution of the Poisson equation we use an integral representation of the electrostatic potential by convolution of the free-space Green’s function with the charge distribution. The convolution equation is solved with the help of the Fast Fourier Transform (FFT). The same algorithm for solution of the 3D Poisson equation is used, for example, in ASTRA

apply(p_array, zstep)

the method is called on every step.

Parameters
  • p_array

  • dz

Returns

el_field(X, Q, gamma, nxyz)
potential(q, steps)
prepare(lat)

method is called at the moment of Physics Process addition to Navigator class.

Parameters

lat

Returns

sym_kernel(ijk2, hxyz)
class ocelot.SpontanRadEffects(K=0.0, lperiod=0.0, type='planar')

Bases: ocelot.cpbd.physics_proc.PhysProc

Effects of the spontaneous radiation: energy loss and quantum diffusion

apply(p_array, dz)

the method is called on every step.

Parameters
  • p_array

  • dz

Returns

energy_loss_und(energy, dz)
sigma_gamma_quant(energy, dz)

rate of energy diffusion

Parameters
  • energy – electron beam energy

  • Kx – undulator parameter

  • lperiod – undulator period

  • dz – length of the

Returns

sigma_gamma/gamma

class ocelot.TDCavity(l=0.0, freq=0.0, phi=0.0, v=0.0, tilt=0.0, eid=None)

Bases: ocelot.cpbd.elements.Element

Transverse deflecting cavity - by default kick in horizontal plane

l - length [m] v - voltage [GV/m] freq - frequency [Hz] phi - phase in [deg] tilt - tilt of cavity in [rad]

class ocelot.TransferMap

Bases: object

TransferMap is a basic linear transfer map for all elements.

apply(prcl_series)
Parameters

prcl_series – can be list of Particles [Particle_1, Particle_2, … ] or ParticleArray

Returns

None

map_x_twiss(tws0)
mul_p_array(rparticles, energy=0.0)
class ocelot.Twiss(beam=None)

Bases: object

class - container for twiss parameters

class ocelot.Undulator(lperiod=0.0, nperiods=0, Kx=0.0, Ky=0.0, field_file=None, eid=None)

Bases: ocelot.cpbd.elements.Element

lperiod - undulator period in [m];

nperiod - number of periods;

Kx - undulator paramenter for vertical field;

Ky - undulator parameter for horizantal field;

field_file - absolute path to magnetic field data;

mag_field - None by default, the magnetic field map function - (Bx, By, Bz) = f(x, y, z) eid - id of undulator.

validate()
class ocelot.UndulatorTestTM(lperiod, Kx, ax=0, ndiv=10)

Bases: ocelot.cpbd.optics.TransferMap

map4undulator(u, z, lperiod, Kx, ax, energy, ndiv)
class ocelot.UnknownElement(l=0, kick=0, xsize=0, ysize=0, volt=0, lag=0, harmon=0, refer=0, vkick=0, hkick=0, eid=None)

Bases: ocelot.cpbd.elements.Element

l - length of lens in [m]

class ocelot.Vcor(l=0.0, angle=0.0, eid=None)

Bases: ocelot.cpbd.elements.RBend

horizontal corrector, l - length of magnet in [m], angle - angle of bend in [rad],

class ocelot.Wake(step=1)

Bases: ocelot.cpbd.physics_proc.PhysProc

The wake field impact on the beam is included as series of kicks. In order to take into account the impact of the wake field on the beam the longitudinal wake function of point charge through the second order Taylor expansion is used. In general case it uses 13 one-dimensional functions to represent the longitudinal component of the wake function for arbitrary sets of the source and the wittness particles near to the reference axis.

w_sampling = 500 - defines the number of the equidistant sampling points for the one-dimensional

wake coefficients in the Taylor expansion of the 3D wake function.

filter_order = 20 - smoothing filter order wake_table = None - wake table [WakeTable()] factor = 1. - scaling coefficient TH - list from WakeTable, (T, H): T- table of wakes coefs, H - matrix of the coefs place in T

add_total_wake(X, Y, Z, q, TH, Ns, NF)
add_wake(I, T)

[x, W] = AddWake(I, T) :param I: wake table in V/C, W in V (R, L, Cinv, nm, W0, N0, W1, N1) :param T: wake table in V/C, W in V :return:

apply(p_array, dz)

the method is called on every step.

Parameters
  • p_array

  • dz

Returns

convolution(xu, u, xw, w)
get_long_wake(current_profile)

method to extract a longitudinal wake from the Table for specific current profile

Parameters

current_profile – 2D array with shape (n, 2) where first column is position and second is a beam current

Returns

wake

prepare(lat)

method is called at the moment of Physics Process addition to Navigator class.

Parameters

lat

Returns

wake_convolution(xb, bunch, xw, wake)
class ocelot.WakeKick(factor=1)

Bases: ocelot.cpbd.wake3D.Wake

apply(p_array, dz)

the method is called on every step.

Parameters
  • p_array

  • dz

Returns

class ocelot.WakeTable(wake_file=None)

Bases: object

WakeTable(wake_file) - load and prepare wake table wake_file - path to the wake table

load_table(wake_file)
process_wake_table(wake_table)
Parameters

wake_file – file name

Returns

(T, H): T- table of wakes coefs, H - matrix of the coefs place in T

read_file(wake_file)
class ocelot.WakeTableDechirperOffAxis(b=0.0005, a=0.01, width=0.02, t=0.00025, p=0.0005, length=1, sigma=3e-05, orient='horz')

Bases: ocelot.cpbd.wake3D.WakeTable

WakeTableDechirperOffAxis() - creates two wake tables for horizontal and vertical corrugated plates. Based on https://doi.org/10.1016/j.nima.2016.09.001 and SLAC-PUB-16881

Parameters
  • b – distance from the plate in [m]

  • a – half gap between plates in [m]

  • width – width of the corrugated structure in [m]

  • t – longitudinal gap in [m]

  • p – period of corrugation in [m]

  • length – length of the corrugated structure in [m]

  • sigma – characteristic (rms) longitudinal beam size in [m]

  • orient – “horz” or “vert” plate orientation

Returns

hor_wake_table, vert_wake_table

calculate_wake_tables(b, a, width, t, p, length, sigma)

Function creates two wake tables for horizontal and vertical corrugated plates

Parameters
  • b – distance from the plate in [m]

  • a – half gap between plates in [m]

  • width – width of the corrugated structure in [m]

  • t – longitudinal gap in [m]

  • p – period of corrugation in [m]

  • length – length of the corrugated structure in [m]

  • sigma – characteristic longitudinal beam size in [m]]

  • filename – save to files if filename is not None

Returns

hor_wake_table, vert_wake_table

ocelot.compensate_chromaticity(lattice, ksi_x_comp=0, ksi_y_comp=0, nsuperperiod=1)

old chromaticity compensation with 2 sextupole families

ocelot.contour_da(track_list, nturns, lvl=0.9)

the function defines contour of DA. If particle “lived” > lvl*nturns then we set up nturns if particle “lived” < lvl*nturns then we set up 0

ocelot.create_track_list(x_array, y_array, p_array, energy=0.0)

the function create list of Pxy

ocelot.ellipse_from_twiss(emit, beta, alpha)
ocelot.fodo_parameters(betaXmean=36.0, L=10.0, verbose=False)
ocelot.freq_analysis(track_list, lat, nturns, harm=True, diap=0.1, nearest=False, nsuperperiods=1)
ocelot.gauss_from_twiss(emit, beta, alpha)
ocelot.generate_parray(sigma_x=0.0001, sigma_px=2e-05, sigma_y=None, sigma_py=None, sigma_tau=0.001, sigma_p=0.0001, chirp=0.01, charge=5e-09, nparticles=200000, energy=0.13, tau_trunc=None, tws=None)

Method to generate ParticleArray with gaussian distribution.

Note: in ParticleArray, {x, px} and {y, py} are canonical coordinates. {tau, p} is not, to make it canonical

the sign of “tau” should be flipped.

Parameters
  • sigma_x – std(x), x is horizontal cartesian coordinate.

  • sigma_px – std(px), ‘px’ is conjugate momentum canonical momentum px/p0.

  • sigma_y – std(y), y is vertical cartesian coordinate.

  • sigma_py – std(py), ‘py’ is canonical momentum py/p0.

  • sigma_tau – std(tau), “tau” = c*t

  • sigma_p – std(p), ‘p’ is canonical momentum E/(c*p0)

  • chirp – energy chirp [unitless], linear correlation - p_i += chirp * tau_i/sigma_tau

  • charge – beam charge in [C], 5e-9 by default

  • nparticles – namber of particles, 200k by default

  • energy – beam energy in [GeV], 0.13 [GeV]

  • tau_trunc – None, if not [float] - truncated gauss distribution in “tau” direction.

  • tws – None, if Twiss obj - the beam is matched to twiss params.

Returns

ParticleArray

ocelot.get_current(p_array, num_bins=200, **kwargs)

Function calculates beam current from particleArray.

Parameters
  • p_array – particleArray

  • charge

    • None, OBSOLETE, charge of the one macro-particle.

    If None, charge of the first macro-particle is used

  • num_bins – number of bins

:return s, I - (np.array, np.array) - beam positions [m] and currents in [A]

ocelot.get_envelope(p_array, tws_i=<ocelot.cpbd.beam.Twiss object>, bounds=None)

Function to calculate twiss parameters form the ParticleArray

Parameters
  • p_array – ParticleArray

  • tws_i – optional, design Twiss,

  • bounds – optional, [left_bound, right_bound] - bounds in units of std(p_array.tau())

Returns

Twiss()

ocelot.get_map(lattice, dz, navi)
ocelot.global_slice_analysis(parray, nparts_in_slice=5000, smooth_param=0.01, filter_base=2, filter_iter=2)

Function to calculate slice parameters

Parameters
  • parray – ParticleArray

  • nparts_in_slice – 5000, nparticles in the slice

  • smooth_param – 0.01, smoothing parameters to calculate the beam current: smooth_param = m_std * np.std(p_array.tau())

  • filter_base – 2, filter parameter in the func: simple_filter

  • filter_iter – 2, filter parameter in the func: simple_filter

Returns

SliceParameters,

ocelot.lattice_transfer_map(lattice, energy)

Function calculates transfer maps, the first and second orders (R, T), for the whole lattice. Second order matrices are attached to lattice object: lattice.T_sym - symmetric second order matrix lattice.T - second order matrix lattice.R - linear R matrix

Parameters
  • lattice – MagneticLattice

  • energy – the initial electron beam energy [GeV]

Returns

R - matrix

ocelot.load_particle_array(filename, print_params=False)

Universal function to load beam file, *.ast (ASTRA), *.fmt1 (CSRTrack) or *.npz format

Note that downloading ParticleArray from the astra file (.ast) and saving it back does not give the same distribution. The difference arises because the array of particles does not have a reference particle, and in this case the first particle is used as a reference.

Parameters

filename – path to file, filename.ast or filename.npz

Returns

ParticleArray

ocelot.match(lat, constr, vars, tw, verbose=True, max_iter=1000, method='simplex', weights=<function weights_default>, vary_bend_angle=False, min_i5=False)

Function to match twiss parameters

Parameters
  • lat – MagneticLattice

  • constr

    dictionary, constrains. Example: ‘periodic’:True - means the “match” function tries to find periodic solution at the ends of lattice:

    constr = {elem1:{‘beta_x’:15, ‘beta_y’:2}, ‘periodic’:True}

    ”hard” constrains on the end of elements (elem1, elem2):

    constr = {elem1:{‘alpha_x’:5, ‘beta_y’:5}, elem2:{‘Dx’:0 ‘Dyp’:0, ‘alpha_x’:5, ‘beta_y’:5}}

    or mixture of “soft” and hard constrains:

    constr = {elem1:{‘alpha_x’:[“>”, 5], ‘beta_y’:5}, elem2:{‘Dx’:0 ‘Dyp’:0, ‘alpha_x’:5, ‘beta_y’:[“>”, 5]}}

    in case one needs global control on beta function, the constrains can be written following way.

    constr = {elem1:{‘alpha_x’:5, ‘beta_y’:5}, ‘global’: {‘beta_x’: [‘>’, 10]}}

  • vars – list of elements e.g. vars = [QF, QD] or it can be initial twiss parameters: vars = [[tws0, ‘beta_x’], [tws0, ‘beta_y’], [tws0, ‘alpha_x’], [tws0, ‘alpha_y’]]

  • tw – initial Twiss

  • verbose – allow print output of minimization procedure

  • max_iter

  • method – string, available ‘simplex’, ‘cg’, ‘bfgs’

  • weights

    function returns weights, for example def weights_default(val):

    if val == ‘periodic’: return 10000001.0 if val == ‘total_len’: return 10000001.0 if val in [‘Dx’, ‘Dy’, ‘Dxp’, ‘Dyp’]: return 10000002.0 if val in [‘alpha_x’, ‘alpha_y’]: return 100007.0 if val in [‘mux’, ‘muy’]: return 10000006.0 if val in [‘beta_x’, ‘beta_y’]: return 100007.0 return 0.0001

  • vary_bend_angle – False, allow to vary “angle” of the dipoles instead of the focusing strength “k1”

  • min_i5 – minimization of the radiation integral I5. Can be useful for storage rings.

Returns

result

ocelot.match_tunes(lat, tw0, quads, nu_x, nu_y, ncells=1, print_proc=0)
ocelot.merger(lat, remaining_types=[], remaining_elems=[], init_energy=0.0)

Function to compress the lattice excluding elements by type or by individual elements

Parameters
  • lat – MagneticLattice

  • remaining_types – list, the type of the elements which needed to be untouched others will be “compress” to Matrix element e.g. [Monitor, Quadrupole, Bend, Hcor]

  • remaining_elems – list of elements (ids (str) or object)

  • init_energy – initial energy

Returns

New MagneticLattice

ocelot.nearest_particle(track_list, xi, yi)
ocelot.save_particle_array(filename, p_array)

Universal function to save beam file, *.ast (ASTRA), *.fmt1 (CSRTrack) or *.npz format

Note that downloading ParticleArray from the astra file (.ast) and saving it back does not give the same distribution. The difference arises because the array of particles does not have a reference particle, and in this case the first particle is used as a reference.

Parameters

filename – path to file, filename.ast or filename.npz

Returns

ParticleArray

ocelot.spectrum(data1D)

input: 1D sample data output: frequency and fourier transform

ocelot.stable_particles(track_list, nturns)
ocelot.track(lattice, p_array, navi, print_progress=True, calc_tws=True, bounds=None)

tracking through the lattice

Parameters
  • lattice – Magnetic Lattice

  • p_array – ParticleArray

  • navi – Navigator

  • print_progress – True, print tracking progress

  • calc_tws – True, during the tracking twiss parameters are calculated from the beam distribution

  • bounds – None, optional, [left_bound, right_bound] - bounds in units of std(p_array.tau())

Returns

twiss_list, ParticleArray. In case calc_tws=False, twiss_list is list of empty Twiss classes.

ocelot.track_nturns(lat, nturns, track_list, nsuperperiods=1, save_track=True, print_progress=True)
ocelot.track_nturns_mpi(mpi_comm, lat, nturns, track_list, errors=None, nsuperperiods=1, save_track=True)
ocelot.tracking_step(lat, particle_list, dz, navi)

tracking for a fixed step dz :param lat: Magnetic Lattice :param particle_list: ParticleArray or Particle list :param dz: step in [m] :param navi: Navigator :return: None

ocelot.twiss(lattice, tws0=None, nPoints=None)

twiss parameters calculation

Parameters
  • lattice – lattice, MagneticLattice() object

  • tws0 – initial twiss parameters, Twiss() object. If None, try to find periodic solution.

  • nPoints – number of points per cell. If None, then twiss parameters are calculated at the end of each element.

Returns

list of Twiss() objects

ocelot.write_lattice(lattice, tws0=None, file_name='lattice.py', remove_rep_drifts=True, power_supply=False)

saves lattice as python imput file lattice - MagneticLattice file_name - name of the file remove_rep_drifts - if True, remove the drifts with the same lengths from the lattice drifts definition