Charged Particle Beam Dynamics (CPBD) module¶
Overview¶
Charged Particle Beam Dynamics module provides features for charged particle (electron) beam optics, including calculating and matching Twiss parameters, single-particle tracking as well as tracking with collective effects (CSR, space charge and wakefields)
Getting started¶
Import OCELOT
from ocelot import *
Define a magnetic lattice
q1 = Quadrupole(l = 0.3, k1 = 5)
q2 = Quadrupole(l = 0.3, k1 = -5)
d = Drift(l = 0.5)
lat = MagneticLattice( (d, q1, d, q2, d, q1, d,q2,d) )
Use twiss()
to find linear optics (Twiss) functions for given initial values
tw0 = Twiss()
tw0.beta_x = 5.
tw0.alpha_x = -0.87
tw0.beta_y = 2.1
tw0.alpha_y = 0.96
tws = twiss(lat, tw0)
Find periodic Twiss solution
tws = twiss(lat)
Find periodic Twiss solution with given longitudinal resolution (500 points)
tws = twiss(lat, nPoints=500)
Plot Twiss parameters
from pylab import *
plot([t.s for t in tws], [t.beta_x for t in tws])
plot([t.s for t in tws], [t.beta_y for t in tws])
Plot Twiss parameters in the lattice display
from ocelot.gui.accelerator import *
plot_opt_func(lat, tws)
show()
Matching¶
-
match
(lattice, constarints, variables[, start=0])¶ lattice a
MagneticLattice
object
Tracking¶
Elements¶
-
class
MagneticLattice
¶
-
class
Drift
¶
-
class
Quadrupole
¶
-
class
Bend
¶ -
same as SBend
-
class
SBend
¶
-
class
RBend
¶
Transfer maps¶
Transfer maps define how the element map acts in tracking. The default transfer map attachment scheme is as follows:
Drifts, Quadrupoles, and bends have first order transfer maps
Sextupoles have a drift-kick-drift map
API documentation¶
-
ocelot.cpbd.optics.
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.cpbd.match.
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