The methods implemented in this module are internal to the ‘lattice’ package.
This is necessary to ensure that the ‘lattice’ package is independent of other
AT packages.
Other Lattice methods are implemented in other AT packages and are available
as soon as the package is imported. The ‘tracking’ and ‘physics’ packages are
automatically imported.
iter – function called as iter(params,*args). It must
return an iterable of Element objects for building
the lattice. It must also fill the params dictionary providing
the Lattice attributes.
params – dictionary of lattice parameters. A custom iterator may
add, remove or modify parameters. Finally, the remaining
parameters will be set as Lattice attributes.
Keyword Arguments:
name='' – Name of the lattice
energy – Energy of the lattice
periodicity=1 – Number of periods
particle='relativistic' – circulating particle. May be
‘relativistic’, ‘electron’, ‘positron’, ‘proton’
or a Particle object
iterator=None – custom iterator (see below)
* – all other keywords will be set as attributes of
the Lattice object
beam_current – Total current in the beam, used for collective
effects
An iterator it is called as it(params,*args) where args
and params are the arguments of the Lattice constructor. It
must yield the AT Elements for building the lattice. It must
also fill its params dictionary argument, which will be used to
set the Lattice attributes.
An iterator can be:
a “generator” which yields elements from scratch.
Examples: a list, or a file iterator,
a “filter” which runs through an input iterator, processes each
element, possibly adds parameters to the params dictionary
and yields the processed elements.
Note
To reduce the inter-package dependencies, some methods of the
lattice object are defined in other AT packages, in the module where
the underlying function is implemented.
Note
It is possible to define a filling pattern for the beam using the
function ring.set_fillingpattern(). The default configuration
(no arguments) is for single bunch and is the one loaded at lattice
initialization. See function help for details.
Changing Lattice.harmonic_number will reset the filling pattern
to its default configuration.
The beam current can be changed with
Lattice.beam_current=current
The filling pattern and beam current are used by collective effects
passmethods.
an integer in the range [-len(ring), len(ring)-1]
selecting the element according to python indexing rules.
As a special case, len(ring) is allowed and refers to the end
of the last element,
an ordered list of such integers without duplicates,
a numpy array of booleans of maximum length len(ring)+1,
where selected elements are True.
Develop a periodical lattice by repeating its elements
self.periodicity times
Parameters:
copy_elements (bool) – If True, deep copies of the elements are used for
the repetition. Otherwise, the original elements are repeated in the
developed lattice.
If False, the modification is done in-place,
If True, return a shallow copy of the lattice. Only the
radiating elements are copied with PassMethod modified.
Caution
a shallow copy means that all non-modified elements are shared
with the original lattice. Any further modification will
affect both lattices.
Turn longitudinal motion on. By default, turn on
radiation in dipoles and quadrupoles, turn on RF cavities, activates
collective effects and other elements acting on momentum.
Modify the lattice in-place or creates a new lattice, depending on the
copy keyword argument.
Syntax using positional arguments:
Parameters:
elem_class –
LongtMotion subclass.
Longitudinal motion is turned on for elements which are
instances of any given elem_class.
In adition to single element classes, a few grouping classes are
available:
LongtMotion: all elements possibly acting on
momentum,
If False, the modification is done in-place,
If True, return a shallow copy of the lattice. Only the
radiating elements are copied with PassMethod modified.
Caution
a shallow copy means that all non-modified elements are shared
with the original lattice. Any further modification will
affect both lattices.
orbit (ndarray | None) – Avoids looking for initial the closed orbit if it is
already known ((6,) array).
keep_lattice (bool) – Assume no lattice change since the previous tracking.
Default: False
Keyword Arguments:
full (bool) – If True, matrices are full 1-turn matrices
at the entrance of each element indexed by refpts. If False
(Default), matrices are between the entrance of the first element and
the entrance of the selected element
Finds the closed orbit in the 4-d transverse phase space by numerically
solving for a fixed point of the one turn map M calculated with
internal_lpass().
under the CONSTANT MOMENTUM constraint dp and with NO
constraint on the 6-th coordinate \(c\tau\)
Important
find_orbit4() imposes a constraint on dp and relaxes the
constraint on the revolution frequency. A physical storage ring does
exactly the opposite: the momentum deviation of a particle on the
closed orbit settles at the value such that the revolution is
synchronous with the RF cavity: \(f_{RF} = h*f_{rev}\)
To impose this artificial constraint in find_orbit4(),
the PassMethod used for any element SHOULD NOT:
Change the longitudinal momentum dp (cavities ,
magnets with radiation)
Have any time dependence (localized impedance, fast kickers etc)
dct (float | None) – Path lengthening. If specified, dp is ignored and
the off-momentum is deduced from the path lengthening.
df (float | None) – Deviation from the nominal RF frequency. If specified,
dp is ignored and the off-momentum is deduced from the frequency
deviation.
orbit (ndarray | None) – Avoids looking for initial the closed orbit if it is
already known ((6,) array). find_orbit4() propagates it to
the specified refpts.
keep_lattice (bool) – Assume no lattice change since the previous tracking.
Default: False
Keyword Arguments:
guess (Orbit) – (6,) initial value for the
closed orbit. It may help convergence. Default: (0, 0, 0, 0, 0, 0)
The requency in the RF cavities \(f_{RF}\) (may be different
from \(f_{RF_0}\)) imposes the synchronous condition
\(c\tau_2 - c\tau_1 = c.h (1/f_{RF} - 1/f_{RF_0})\),
The algorithm numerically calculates the
6-by-6 Jacobian matrix J6. In order for (J-E) matrix
to be non-singular it is NECESSARY to use a realistic
PassMethod for cavities with non-zero momentum kick
(such as RFCavityPass).
find_orbit6() can find orbits with radiation.
In order for the solution to exist the cavity must supply an
adequate energy compensation.
In the simplest case of a single cavity, it must have
Voltage set so that \(V > E_{loss}\), the energy loss
per turn
There is a family of solutions that correspond to different RF
buckets They differ in the 6-th coordinate by
\(c*Nb/f_{RF}\), Nb = 0:h-1
The value of the 6-th coordinate found at the cavity gives
the equilibrium RF phase. If there is no radiation it is 0.
dp, dct and df arguments are applied with respect
to the NOMINAL on-momentum frequency. They overwrite
exisiting frequency offsets
orbit (ndarray) – Avoids looking for initial the closed orbit if it is
already known ((6,) array). find_sync_orbit() propagates it
to the specified refpts.
keep_lattice (bool) – Assume no lattice change since the previous tracking.
Default: False
Keyword Arguments:
guess (Orbit) – (6,) initial value for the
closed orbit. It may help convergence. Default: (0, 0, 0, 0, 0, 0)
Finds the closed orbit, synchronous with the RF cavity (first 5
components of the phase space vector) by numerically solving for a fixed
point of the one turn map M calculated with internal_lpass()
under the constraint \(dc\tau = c\tau_2 - c\tau_1 = c/f_{rev} -
c/f_{rev_0}\), where \(f_{rev_0} = f_{RF_0}/h\) is the design revolution
frequency and \(f_{rev} = (f_{RF_0} + df_{RF})/h\) is the imposed
revolution frequency.
Important
find_sync_orbit() imposes a constraint \(c\tau_2 -
c\tau_1\) and \(dp_2 = dp_1\) but no constraint on the value of
\(dp_2\) or \(dp_1\). The algorithm assumes time-independent
fixed-momentum ring to reduce the dimensionality of the problem.
To impose this artificial constraint in find_sync_orbit(),
the PassMethod used for any element SHOULD NOT:
Change the longitudinal momentum dp (cavities ,
magnets with radiation)
Have any time dependence (localized impedance, fast kickers etc)
dp (float | None) – Momentum deviation. Defaults to None
df (float | None) – Deviation from the nominal RF frequency. If specified,
dct is ignored and the off-momentum is deduced from the frequency
deviation.
orbit (ndarray | None) – Avoids looking for initial the closed orbit if it is
already known ((6,) array). find_sync_orbit() propagates it
to the specified refpts.
keep_lattice (bool) – Assume no lattice change since the previous tracking.
Default: False
Keyword Arguments:
guess (Orbit) – (6,) initial value for the
closed orbit. It may help convergence. Default: (0, 0, 0, 0, 0, 0)
verbose (bool | None) – Print out some information
divider (int | None) – Value of the divider used in
GridMode.RECURSIVE boundary search
shift_zero (float | None) – Epsilon offset applied on all 6 coordinates
start_method (str | None) – Python multiprocessing start method. The default
None uses the python default that is considered safe.
Available parameters: 'fork', 'spawn', 'forkserver'.
The default for linux is 'fork', the default for macOS and
Windows is 'spawn'. 'fork' may be used on macOS to speed up
the calculation or to solve runtime errors, however it is
considered unsafe.
Returns:
boundary – (2,n) array: 2D acceptance
survived – (2,n) array: Coordinates of surviving particles
tracked – (2,n) array: Coordinates of tracked particles
In case of multiple refpts, return values are lists of arrays, with one
array per ref. point.
When``use_mp=True`` all the available CPUs will be used.
This behavior can be changed by setting
at.DConstant.patpass_poolsize to the desired value
When multiple refpts are provided particles are first
projected to the beginning of the ring with tracking. Then,
all particles are tracked up to nturns. This allows to
do most of the work in a single function call and allows for
full parallelization.
Solves the Haissinski formula and returns the bunch length and energy
spread for given bunch current and \(Z/n\). If both zn and
bunch_curr are None, zero current case, otherwise both are needed
for the calculation
Parameters:
ring – ring use for tracking
Keyword Arguments:
zn=None – \(Z/n\) for the full ring
bunch_curr=None – Bunch current
espread=None – Energy spread, if None use lattice parameter
Returns a bool array of element indices, selecting ring elements.
Deprecated: get_cells(ring,refpts) is
ring.bool_refpts(refpts) except for str arguments:
get_cells(ring,attrname[,attrvalue]) is
ring.bool_refpts(checkattr(strkey[,attrvalue]))
use_mp – Use python multiprocessing (patpass(),
default use lattice_pass()). In case multiprocessing
is not enabled, grid_mode is forced to
GridMode.RECURSIVE (most efficient in single core)
verbose – Print out some information
divider – Value of the divider used in
GridMode.RECURSIVE boundary search
shift_zero – Epsilon offset applied on all 6 coordinates
start_method – Python multiprocessing start method. The default
None uses the python default that is considered safe.
Available parameters: 'fork', 'spawn', 'forkserver'.
The default for linux is 'fork', the default for macOS and
Windows is 'spawn'. 'fork' may be used on macOS to speed up
the calculation or to solve runtime errors, however it is considered
unsafe.
Returns:
boundary – (len(refpts),2) array: 1D acceptance
survived – (n,) array: Coordinates of surviving particles
tracked – (n,) array: Coordinates of tracked particles
In case of multiple tracked and survived are lists of arrays,
with one array per ref. point.
Note
When``use_mp=True`` all the available CPUs will be used.
This behavior can be changed by setting
at.DConstant.patpass_poolsize to the desired value
When multiple refpts are provided particles are first
projected to the beginning of the ring with tracking. Then,
all particles are tracked up to nturns. This allows to
do most of the work in a single function call and allows for
full parallelization.
Computes the touschek lifetime using the Piwinski formula
Parameters:
ring – ring use for tracking
emity – verticla emittance
bunch_curr – bunch current
Keyword Arguments:
emitx=None – horizontal emittance
sigs=None – rms bunch length
sigp=None – energy spread
zn=None – full ring \(Z/n\)
momap=None – momentum aperture, has to be consistent with
refpts if provided the momentum aperture is
not calculated
refpts=None – refpts where the momentum aperture is calculated,
the default is to compute it for all elements in the
ring, len(refpts)>2 is required
resolution – minimum distance between 2 grid points, default=1.0e-3
amplitude – max. amplitude for RADIAL and CARTESIAN or
initial step in RECURSIVE
default = 0.1
nturns=1024 – Number of turns for the tracking
dp=None – static momentum offset
offset – initial orbit. Default: closed orbit
grid_mode – at.GridMode.CARTESIAN/RADIAL track full vector
(default). at.GridMode.RECURSIVE: recursive search
use_mp=False – Use python multiprocessing (patpass, default use
lattice_pass). In case multi-processing is not
enabled GridMode is forced to
RECURSIVE (most efficient in single core)
verbose=False – Print out some inform
epsrel (epsabs,) – integral absolute and relative tolerances
Returns:
tl – touschek lifetime, double expressed in seconds
ma – momentum aperture (len(refpts), 2) array
refpts – refpts used for momentum aperture calculation
(len(refpts), ) array
Note
When``use_mp=True`` all the available CPUs will be used.
This behavior can be changed by setting
at.DConstant.patpass_poolsize to the desired value
When multiple refpts are provided particles are first
projected to the beginning of the ring with tracking. Then,
all particles are tracked up to nturns. This allows to
do most of the work in a single function call and allows for
full parallelization.
ring (Lattice) – Lattice description (ring.is_6d must be
False)
dp (float | None) – Momentum deviation. Defaults to None
keep_lattice (bool) – Assume no lattice change since the previous tracking.
Default: False
fit_order (int | None) – Maximum momentum compaction factor order to be fitted.
Default to 1, corresponding to the first-order momentum compaction factor.
n_step (int | None) – Number of different calculated momentum deviations to be fitted
with a polynomial.
Default to 2.
use_mp – Use python multiprocessing (patpass(),
default use lattice_pass()). In case multiprocessing is
not enabled, grid_mode is forced to
GridMode.RECURSIVE (most efficient in single core)
verbose – Print out some information
divider – Value of the divider used in
GridMode.RECURSIVE boundary search
shift_zero – Epsilon offset applied on all 6 coordinates
start_method – Python multiprocessing start method. The default
None uses the python default that is considered safe.
Available parameters: 'fork', 'spawn', 'forkserver'.
The default for linux is 'fork', the default for macOS and
Windows is 'spawn'. 'fork' may be used on macOS to speed up
the calculation or to solve runtime errors, however it is considered
unsafe.
Returns:
boundary – (len(refpts),2) array: 1D acceptance
survived – (n,) array: Coordinates of surviving particles
tracked – (n,) array: Coordinates of tracked particles
In case of multiple tracked and survived are lists of arrays,
with one array per ref. point.
Note
When``use_mp=True`` all the available CPUs will be used.
This behavior can be changed by setting
at.DConstant.patpass_poolsize to the desired value
When multiple refpts are provided particles are first
projected to the beginning of the ring with tracking. Then,
all particles are tracked up to nturns. This allows to
do most of the work in a single function call and allows for
full parallelization.
an integer in the range [-len(ring), len(ring)-1]
selecting the element according to python indexing rules.
As a special case, len(ring) is allowed and refers to the end
of the last element,
an ordered list of such integers without duplicates,
a numpy array of booleans of maximum length len(ring)+1,
where selected elements are True.
Computes the touschek scattering using the Piwinski formula
Parameters:
ring – ring use for tracking
emity – verticla emittance
bunch_curr – bunch current
Keyword Arguments:
emitx=None – horizontal emittance
sigs=None – rms bunch length
sigp=None – energy spread
zn=None – full ring \(Z/n\)
momap=None – momentum aperture, has to be consistent with
refpts if provided the momentum aperture
is not calculated
refpts=None – refpts where the momentum aperture is calculated,
the default is to compute it for all elements in the
ring, len(refpts)>2 is required
resolution – minimum distance between 2 grid points, default=1.0e-3
amplitude – max. amplitude for RADIAL and CARTESIAN or
initial step in RECURSIVE
default = 0.1
nturns=1024 – Number of turns for the tracking
dp=None – static momentum offset
grid_mode – at.GridMode.CARTESIAN/RADIAL track full vector
(default). at.GridMode.RECURSIVE: recursive search
use_mp=False – Use python multiprocessing (patpass, default use
lattice_pass). In case multi-processing is not
enabled GridMode is forced to
RECURSIVE (most efficient in single core)
verbose=False – Print out some inform
epsrel (epsabs,) – integral absolute and relative tolerances
use_mp – Use python multiprocessing (patpass(),
default use lattice_pass()). In case multiprocessing
is not enabled, grid_mode is forced to
GridMode.RECURSIVE (most efficient in single core)
verbose – Print out some information
divider – Value of the divider used in
GridMode.RECURSIVE boundary search
shift_zero – Epsilon offset applied on all 6 coordinates
start_method – Python multiprocessing start method. The default
None uses the python default that is considered safe.
Available parameters: 'fork', 'spawn', 'forkserver'.
The default for linux is 'fork', the default for macOS and
Windows is 'spawn'. 'fork' may be used on macOS to speed up
the calculation or to solve runtime errors, however it is considered
unsafe.
Returns:
boundary – (len(refpts),2) array: 1D acceptance
survived – (n,) array: Coordinates of surviving particles
tracked – (n,) array: Coordinates of tracked particles
In case of multiple tracked and survived are lists of arrays,
with one array per ref. point.
Note
When``use_mp=True`` all the available CPUs will be used.
This behavior can be changed by setting
at.DConstant.patpass_poolsize to the desired value
When multiple refpts are provided particles are first
projected to the beginning of the ring with tracking. Then,
all particles are tracked up to nturns. This allows to
do most of the work in a single function call and allows for
full parallelization.
lattice_pass() tracks particles through each element of a lattice
calling the element-specific tracking function specified in the Element’s
PassMethod field.
r_in – (6, N) array: input coordinates of N particles.
r_in is modified in-place and reports the coordinates at
the end of the element. For the best efficiency, r_in
should be given as F_CONTIGUOUS numpy array.
keep_lattice (bool) – Use elements persisted from a previous
call. If True, assume that the lattice has not changed
since the previous call.
keep_counter (bool) – Keep the turn number from the previous
call.
turn (int) – Starting turn number. Ignored if
keep_counter is True. The turn number is necessary to
compute the absolute path length used in RFCavityPass.
losses (bool) – Boolean to activate loss maps output
omp_num_threads (int) – Number of OpenMP threads
(default: automatic)
The following keyword arguments overload the Lattice values
energy (Optiona[float]) – lattice energy. Default 0.
unfold_beam (bool) – Internal beam folding activate, this
assumes the input particles are in bucket 0, works only
if all bucket see the same RF Voltage.
Default: True
If energy is not available, relativistic tracking if forced,
rest_energy is ignored.
Returns:
r_out – (6, N, R, T) array containing output coordinates of N particles
at R reference points for T turns.
loss_map – If losses is True: dictionary with the
following key:
islost
(npart,) bool array indicating lost particles
turn
(npart,) int array indicating the turn at
which the particle is lost
element
((npart,) int array indicating the element at
which the particle is lost
coord
(6, npart) float array giving the coordinates at
which the particle is lost (zero for surviving
particles)
Note
lattice_pass(lattice,r_in,refpts=len(line)) is the same as
lattice_pass(lattice,r_in) since the reference point
len(line) is the exit of the last element.
lattice_pass(lattice,r_in,refpts=0) is a copy of r_in
since the reference point 0 is the entrance of the first element.
To resume an interrupted tracking (for instance to get intermediate
results), one must use one of the turn or keep_counter
keywords to ensure the continuity of the turn number.
For multiparticle tracking with large number of turn the size of
r_out may increase excessively. To avoid memory issues
lattice_pass(lattice,r_in,refpts=None) can be used.
An empty list is returned and the tracking results of the last turn
are stored in r_in.
To model buckets with different RF voltage unfold_beam=False
has to be used. The beam can be unfolded using the function
unfold_beam(). This function takes into account
the true voltage in each bucket and distributes the particles in the
bunches defined by ring.fillpattern using a 6D orbit search.
an integer in the range [-len(ring), len(ring)-1]
selecting the element according to python indexing rules.
As a special case, len(ring) is allowed and refers to the end
of the last element,
an ordered list of such integers without duplicates,
a numpy array of booleans of maximum length len(ring)+1,
where selected elements are True.
df (float) – Deviation of RF frequency. Defaults to
None
orbit (Orbit) – Avoids looking for the closed orbit if it is
already known ((6,) array)
get_chrom (bool) – Compute chromaticities. Needs computing
the tune at 2 different momentum deviations around the central one.
get_w (bool) – Computes chromatic amplitude functions
(W, WP) [4], and derivatives of the dispersion and twiss parameters
versus dp. Needs to compute the optics at 2 different momentum
deviations around the central one.
keep_lattice (bool) – Assume no lattice change since the
previous tracking. Defaults to False
an integer in the range [-len(ring), len(ring)-1]
selecting the element according to python indexing rules.
As a special case, len(ring) is allowed and refers to the end
of the last element,
an ordered list of such integers without duplicates,
a numpy array of booleans of maximum length len(ring)+1,
where selected elements are True.
an integer in the range [-len(ring), len(ring)-1]
selecting the element according to python indexing rules.
As a special case, len(ring) is allowed and refers to the end
of the last element,
an ordered list of such integers without duplicates,
a numpy array of booleans of maximum length len(ring)+1,
where selected elements are True.
The file format is indicated by the filepath extension. The file name is stored in
the in_file Lattice attribute. The selected variable, if relevant, is stored
in the use Lattice attribute.
Modify selected elements, in-place or in a lattice copy
Parameters:
elem_modify (Callable) – element selection function.
If elem_modify(elem) returns None, the element is
unchanged. Otherwise, elem_modify(elem) must return a
dictionary of attribute name and values, to be set to elem.
copy (bool | None) – If True, return a shallow copy of the lattice.
Only the modified elements are copied.
If False, the modification is done in-place
Returns:
newring – New lattice if copy is True, None
if copy is False
Plot the resulting longitudinal Hamiltonian of a ring (defining the RF
bucket). The Hamiltonian is calculated by summing all the cavities in the
ring. Harmonic cavities are supported by the function.
A perfectly tuned lattice is assumed, the cavities’ frequency is nominal
and the TimeLag is set in a way ensuring ct=0 for the synchronous phase
by using set_cavity_phase().
Parameters:
ring – Lattice description
ct_range (tuple) – Forced lower and upper ct values for the plot.
Default to \(\pm 1.1 \times C / (2h)\)
dp_range (tuple) – Forced lower and upper dp values for the plot.
Default to twice the RF acceptance of the bucket.
num_points (int) – Number of points for 1D grid (ct/dp)
Default to 400.
num_levels (int) – Number of levels for contour plot. Odd number of
levels allow to center the colorbar around 0. Default to 41.
plot_separatrix (bool) – Flag to plot the separatrix contour
(\(\mathcal{H}=0\)).
use_mp (bool) – Use python multiprocessing (patpass(),
default use lattice_pass()). In case multiprocessing is not
enabled, grid_mode is forced to GridMode.RECURSIVE
(most efficient in single core)
start_method (str) – Python multiprocessing start method. The default
None uses the python default that is considered safe.
Available parameters: ‘fork’, ‘spawn’, ‘forkserver’.
The default for linux is ‘fork’, the default for macOS and
Windows is ‘spawn’. ‘fork’ may be used for macOS to speed up
the calculation or to solve runtime errors, however it is
considered unsafe.
Returns:
boundary – (2,n) array: 2D acceptance
tracked – (2,n) array: Coordinates of tracked particles
survived – (2,n) array: Coordinates of surviving particles
Removes all elements with an IdentityPass PassMethod and merges
compatible consecutive elements.
The “reduced” lattice has the same optics as the original one, but
fewer elements. Compatible elements must have the same PolynomA,
PolynomB and bending radius, so that the optics is preserved. But
magnet misalignments are not preserved, so this method should be
applied to lattices without errors.
Keyword Arguments:
keep (Refpts) – Keep the selected elements, even with
IdentityPass PassMethod. Default: keep Monitor
and RFCavity elements
This method allows to repeat the lattice n times.
If n does not divide ring.periodicity, the new ring
periodicity is set to 1, otherwise it is set to
ring.periodicity /= n.
copy_elements (bool) – If True, deep copies of the lattice are used for
the repetition. Otherwise, the original elements are repeated in the
developed lattice.
Insert elements at selected locations in the lattice
Parameters:
break_s – location or array of locations of breakpoints
break_elems – elements to be inserted at breakpoints (array of
elements as long as break_s or single element
duplicated as necessary). Default: Marker(‘sbreak’)
Returns:
newring – A new lattice with new elements inserted at breakpoints
array (bool | None) – If False, the value is applied as described for
set_rf_voltage, set_rf_timelag and set_rf_frequency
If True, directly apply the value to the selected
cavities. The value must be broadcastable to the number
of cavities.
copy (bool | None) – If True, returns a shallow copy of ring with new
cavity elements. Otherwise, modify ring in-place
Adjust the TimeLag attribute of RF cavities based on frequency,
voltage and energy loss per turn, so that the synchronous phase is zero.
An error occurs if all cavities do not have the same frequency.
Warning
This function changes the time reference, this should be avoided
Function to generate the filling pattern lof the ring.
The filling pattern is computed as:
bunches/numpy.sum(bunches)
This function also generates the bunch spatial distribution
accessible with Lattice.bunch_spos
Keyword Arguments:
bunches – integer or array of positive double or bool to define
the bunch distribution.
For scalar input, equidistant bunches are assumed.
ring.harmonic_number has to be a multiple of
bunches.
For array input the condition
len(bunches)==ring.harmonic_number is required.
(default=1, single bunch configuration).
If False (default), frequency
is applied to the selected cavities with the lowest frequency. The
frequency of all the other selected cavities is scaled by the same
ratio.
If True, directly apply frequency to the selected
cavities. The value must be broadcastable to the number of cavities.
copy (Optional[bool]) – If True, returns a shallow copy of
ring with new cavity elements. Otherwise (default), modify
ring in-place
Scales magnet strengths with local energy to cancel the closed orbit
and optics errors due to synchrotron radiations. The dipole bending angle
is changed by changing the reference momentum. This is the ideal
tapering scheme where magnets and multipoles components (PolynomB and
PolynomA) are scaled individually.
Warning
This method works only for lattices without errors and
corrections: if not all corrections and field errors will also be
scaled !!!
track_function() tracks particles through each element of a
lattice or throught a single Element calling the element-specific
tracking function specified in the Element’s PassMethod field.
r_in – (6, N) array: input coordinates of N particles.
r_in is modified in-place only if in_place is
True and reports the coordinates at
the end of the element. For the best efficiency, r_in
should be given as F_CONTIGUOUS numpy array.
in_place (bool) – If True r_in is modified in-place and
reports the coordinates at the end of the element.
(default: False)
keep_lattice (bool) – Use elements persisted from a previous
call. If True, assume that the lattice has not changed
since the previous call.
keep_counter (bool) – Keep the turn number from the previous
call.
turn (int) – Starting turn number. Ignored if
keep_counter is True. The turn number is necessary to
compute the absolute path length used in RFCavityPass.
losses (bool) – Boolean to activate loss maps output
omp_num_threads (int) – Number of OpenMP threads
(default: automatic)
use_mp (bool) – Flag to activate multiprocessing (default: False)
pool_size – number of processes used when
use_mp is True. If None, min(npart,nproc)
is used. It can be globally set using the variable
at.lattice.DConstant.patpass_poolsize
start_method – python multiprocessing start method.
None uses the python default that is considered safe.
Available values: 'fork', 'spawn', 'forkserver'.
Default for linux is 'fork', default for macOS and Windows is
'spawn'. 'fork' may be used on macOS to speed up the
calculation or to solve Runtime Errors, however it is considered
unsafe. Used only when use_mp is True. It can be globally
set using the variable at.lattice.DConstant.patpass_startmethod
The following keyword arguments overload the lattice values
energy (Optiona[float]) – lattice energy. Default 0.
unfold_beam (bool) – Internal beam folding activate, this
assumes the input particles are in bucket 0, works only
if all bucket see the same RF Voltage.
Default: True
If energy is not available, relativistic tracking if forced,
rest_energy is ignored.
Returns:
r_out – (6, N, R, T) array containing output coordinates of N particles
at R reference points for T turns
trackparam – A dictionary containing tracking input parameters with the
following keys:
npart
number of particles
rout
final particle coordinates
turn
starting turn
refpts
array of index where particle coordinate are saved
(only for lattice tracking)
nturns
number of turn
trackdata – A dictionary containing tracking data with the following
keys:
loss_map:
recarray containing the loss_map (only for lattice
tracking)
The loss_map is filled only if losses is True,
it contains the following keys:
islost
(npart,) bool array indicating lost particles
turn
(npart,) int array indicating the turn at
which the particle is lost
element
(npart,) int array indicating the element at
which the particle is lost
coord
(npart, 6) float array giving the coordinates at
which the particle is lost (zero for surviving
particles)
Note
track_function(lattice,r_in,refpts=len(line)) is the same
as track_function(lattice,r_in) since the reference point
len(line) is the exit of the last element.
track_function(lattice,r_in,refpts=0) is a copy of r_in
since the reference point 0 is the entrance of the first element.
To resume an interrupted tracking (for instance to get intermediate
results), one must use one of the turn or keep_counter
keywords to ensure the continuity of the turn number.
For multiparticle tracking with large number of turn the size of
r_out may increase excessively. To avoid memory issues
track_function(lattice,r_in,refpts=None,in_place=True)
can be used. An empty list is returned and the tracking results of
the last turn are stored in r_in.
To model buckets with different RF voltage unfold_beam=False
has to be used. The beam can be unfolded using the function
unfold_beam(). This function takes into account
the true voltage in each bucket and distributes the particles in the
bunches defined by ring.fillpattern using a 6D orbit search.