Optics calculations
Optics calculations can be performed via the compute
module.
Closed orbit search
Using compute.closed_orbit
we can perform closed orbit search for a given lattice:
[1]:
from importlib import resources
from dipas.build import from_file
from dipas.compute import closed_orbit, linear_closed_orbit
from dipas.elements import Kicker, HKicker, VKicker
import dipas.test.sequences
with resources.path(dipas.test.sequences, 'cryring.seq') as path:
lattice = from_file(path)
thin = lattice.makethin({Kicker: 1})
print(closed_orbit(thin))
tensor([[0.],
[0.],
[0.],
[0.],
[0.],
[0.]])
As can be seen from the above example we first need to convert all elements that don’t support thick tracking (the kicker magnets) to thin elements because the closed orbit search is performed by tracking the closed orbit through the lattice. Since all kickers are turned off (kick == 0
) the closed orbit is just zero. Let’s add some kicks:
[2]:
import random
random.seed(1)
for kicker in lattice[HKicker] + lattice[VKicker]:
kicker.kick = random.uniform(-0.001, 0.001)
thin_with_kicks = lattice.makethin({Kicker: 1}) # Use 'makethin' again to make the added kicks effective.
print(closed_orbit(thin_with_kicks, order=1).flatten())
print(linear_closed_orbit(thin_with_kicks).flatten())
tensor([0.0012, 0.0005, 0.0055, 0.0002, 0.0000, 0.0000])
tensor([0.0012, 0.0005, 0.0055, 0.0002, 0.0000, 0.0000])
One important thing to note is that we need to assign the kicks to the original lattice, since the thin
version doesn’t contain the original kickers anymore. Then for the new kicker values we need to makethin
the lattice
again before performing the closed orbit search. This seems somewhat repetitive but it is important in order to maintain the relation between thick and thin elements. Especially if optimization parameters are involved, it is important to always makethin
the
original lattice (which stores the optimization parameters) in order to always get the up-to-date values of the optimization parameters. linear_closed_orbit
computes the closed orbit directly from first order transfer maps (using a slightly different algorithm).
Twiss computation
We can use compute.twiss
in order to compute lattice functions as well as phase advance and tunes:
[3]:
import dipas.compute as compute
twiss = compute.twiss(thin) # Returns a dict.
print(list(twiss), end='\n\n')
print('Q1:', twiss['Q1'], ', Q2:', twiss['Q2'], end='\n\n')
print('Coupling Matrix:', twiss['coupling_matrix'], sep='\n', end='\n\n')
print('One-Turn Matrix:', twiss['one_turn_matrix'], sep='\n', end='\n\n')
print('Lattice:', twiss['lattice'].columns, sep='\n')
['lattice', 'coupling_matrix', 'Q1', 'Q2', 'one_turn_matrix']
Q1: tensor(2.4200) , Q2: tensor(2.4200)
Coupling Matrix:
tensor([[0., 0.],
[0., 0.]])
One-Turn Matrix:
tensor([[-8.7631e-01, 9.2467e-01, 0.0000e+00, 0.0000e+00, 2.7810e-17,
-8.3255e+00],
[-2.5099e-01, -8.7631e-01, 0.0000e+00, 0.0000e+00, -3.4146e-17,
-1.1137e+00],
[ 0.0000e+00, 0.0000e+00, -8.7631e-01, 1.0994e+00, 0.0000e+00,
0.0000e+00],
[ 0.0000e+00, 0.0000e+00, -2.1111e-01, -8.7631e-01, 0.0000e+00,
0.0000e+00],
[ 1.1137e+00, 8.3255e+00, 0.0000e+00, 0.0000e+00, 1.0000e+00,
3.1825e+02],
[ 1.7416e-17, 1.2832e-17, 0.0000e+00, 0.0000e+00, 1.7203e-33,
1.0000e+00]])
Lattice:
Index(['x', 'px', 'y', 'py', 'bx', 'ax', 'mx', 'by', 'ay', 'my', 'dx', 'dpx',
'dy', 'dpy'],
dtype='object')
Here twiss['lattice']
is a pandas.DataFrame
with element labels as index and lattice functions as columns.
Transfer maps
By using compute.transfer_maps
we can compute the transfer maps along the lattice. The parameter method
let’s us specify how the transfer maps are computed. The following options are available:
method='local'
: Compute the local maps of elements, including closed orbit contribution.method='accumulate'
: Compute the cumulative transfer maps w.r.t. to the start of the lattice.method='reduce'
: Compute the combined transfer map for the whole segment.
[4]:
from dipas.compute import transfer_maps
maps = dict(transfer_maps(thin, method='accumulate', labels=True))
print(*maps[lattice[HKicker, 0].label], sep='\n')
tensor([[0.],
[0.],
[0.],
[0.],
[0.],
[0.]])
tensor([[1.0000, 0.8328, 0.0000, 0.0000, 0.0000, 0.0000],
[0.0000, 1.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[0.0000, 0.0000, 1.0000, 0.8328, 0.0000, 0.0000],
[0.0000, 0.0000, 0.0000, 1.0000, 0.0000, 0.0000],
[0.0000, 0.0000, 0.0000, 0.0000, 1.0000, 6.1278],
[0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 1.0000]])
tensor([[[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, -1.2038],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, -1.2038, 0.0000, 0.0000, 0.0000, 0.0000]],
[[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000]],
[[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, -1.2038],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, -1.2038, 0.0000, 0.0000]],
[[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000]],
[[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, -1.2038, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, -1.2038, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, -26.5733]],
[[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000]]])
Since the kicker represents the full element as a sequence of drift spaces and thin slices, its first order map (i.e. the transfer matrix) is not the identity matrix.
If we want the element local maps instead we can use method='local'
(here 0.3755
is the length of the kicker):
[5]:
maps = dict(transfer_maps(thin, method='local', labels=True))
print('Length of Kicker:', lattice[HKicker, 0].l)
print('Transfer Map:', *maps[lattice[HKicker, 0].label], sep='\n')
Length of Kicker: tensor(0.3755)
Transfer Map:
tensor([[0.],
[0.],
[0.],
[0.],
[0.],
[0.]])
tensor([[1.0000, 0.3755, 0.0000, 0.0000, 0.0000, 0.0000],
[0.0000, 1.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[0.0000, 0.0000, 1.0000, 0.3755, 0.0000, 0.0000],
[0.0000, 0.0000, 0.0000, 1.0000, 0.0000, 0.0000],
[0.0000, 0.0000, 0.0000, 0.0000, 1.0000, 2.7629],
[0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 1.0000]])
tensor([[[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, -0.5428],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, -0.5428, 0.0000, 0.0000, 0.0000, 0.0000]],
[[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000]],
[[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, -0.5428],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, -0.5428, 0.0000, 0.0000]],
[[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000]],
[[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, -0.5428, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, -0.5428, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, -11.9816]],
[[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000],
[ 0.0000, 0.0000, 0.0000, 0.0000, 0.0000, 0.0000]]])
Orbit Response Matrix
Using compute.orm
we can compute the orbit response matrix for a given lattice. We need to specify the kickers and monitors to be used, which can be done in a similar way as for selecting lattice elements in general: either we can specify an identifier that selects multiple elements directly, such as a lattice element type or a regex, or we can specify a list of single element identifiers, such as unambiguous labels for example. Let’s compute the horizontal ORM for one of the example
lattices:
[6]:
from importlib import resources
from dipas.build import from_file
from dipas.compute import orm
from dipas.elements import HKicker, HMonitor
import dipas.test.sequences
with resources.path(dipas.test.sequences, 'cryring.seq') as path:
lattice = from_file(path)
orm_x, orm_y = orm(lattice, kickers=HKicker, monitors=HMonitor)
Here we don’t need to call makethin
beforehand because this will be done inside the orm
function. This is necessary because the orm
function will temporarily vary the kicker strengths and, as explained above, for each change to the original lattice we need to create a new thin version (i.e. changes to the original lattice are not automatically mapped to any thin versions that have been created before).
[7]:
print('ORM.shape: ', orm_x.shape)
print('Number of HKickers and HMonitors: ', (len(lattice[HKicker]), len(lattice[HMonitor])), end='\n\n')
print('ORM-X\n\n', orm_x, end='\n\n')
print('ORM-Y\n\n', orm_y)
ORM.shape: torch.Size([9, 12])
Number of HKickers and HMonitors: (12, 9)
ORM-X
tensor([[ 1.6402, 1.5428, 0.8887, 2.0276, 1.6111, 1.0921, -2.9513, 3.2348,
-2.7059, 1.2091, 2.1282, 2.0422],
[ 0.8140, 1.0098, 1.7728, 0.9366, 1.2233, 1.8082, -1.3954, 0.5891,
0.4745, -1.3684, -0.1676, 0.0053],
[ 0.4583, 0.6885, 1.6933, 1.2629, 0.8969, 1.6400, -0.7518, -0.1673,
1.1597, -1.7372, -0.6953, -0.4921],
[-0.3921, -0.8058, -2.7059, -2.4109, -0.0308, 1.2010, 1.1115, 1.0460,
-2.6066, 3.2348, 1.6820, 1.3166],
[ 1.1925, 1.3288, 1.7242, 0.4483, -1.3899, -2.0740, 1.5688, 1.4899,
-0.4258, -0.7906, 0.5093, 0.6296],
[ 1.4039, 1.7722, 3.2348, 1.8029, -1.3615, -2.7341, 0.9301, 0.8887,
1.0460, -2.6066, -0.4424, -0.1171],
[-1.9135, -2.1046, -2.6066, -0.5499, 2.2668, 3.2892, -2.6391, 1.2091,
0.8887, 1.0460, -0.9556, -1.1243],
[-0.4376, -0.6838, -1.7705, -1.3686, 0.2506, 1.0372, 0.2900, -1.1510,
1.8197, 1.8726, 0.7970, 0.5795],
[ 1.7387, 1.6842, 1.0460, -0.8998, -2.3613, -2.6678, 3.4040, -2.7059,
1.2091, 0.8887, 2.0117, 1.9636]])
ORM-Y
tensor([[0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.],
[0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.],
[0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.],
[0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.],
[0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.],
[0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.],
[0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.],
[0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.],
[0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0., 0.]])