"""Lattice objectThe methods implemented in this module are internal to the 'lattice' package.This is necessary to ensure that the 'lattice' package is independent of otherAT packages.Other Lattice methods are implemented in other AT packages and are availableas soon as the package is imported. The 'tracking' and 'physics' packages areautomatically imported.As an example, see the at.physics.orbit module"""from__future__importannotations__all__=["Lattice","Filter","type_filter","params_filter","lattice_filter","elem_generator","no_filter",]importsysimportcopyimportmathifsys.version_info.minor<9:fromtypingimportCallable,Iterable,GeneratorSupportsIndex=intelse:fromtypingimportSupportsIndexfromcollections.abcimportCallable,Iterable,Generatorimportwarningsimportnumpyasnpfrom.importelementsaseltfrom.elementsimportElementfrom.exceptionsimportAtError,AtWarningfrom.particle_objectimportParticlefrom.utilsimportget_s_pos,get_elements,get_value_refpts,Refpts,set_value_refpts# noinspection PyProtectedMemberfrom.utilsimportget_uint32_index,get_bool_index,_refcount,Uint32Refptsfrom.utilsimportrefpts_iterator,checktype,get_geometryfrom.transformationimportset_rotation,set_tilt,set_shiftfrom..constantsimportclight,e_mass_TWO_PI_ERROR=1.0e-4Filter=Callable[...,Iterable[Element]]_DEFAULT_PASS={False:(("cavity_pass",elt.RFCavity,"auto"),("dipole_pass",elt.Dipole,"auto"),("quadrupole_pass",elt.Quadrupole,"auto"),("wiggler_pass",elt.Wiggler,"auto"),("sextupole_pass",elt.Sextupole,"auto"),("octupole_pass",elt.Octupole,"auto"),("multipole_pass",elt.Multipole,"auto"),("collective_pass",elt.Collective,"auto"),("diffusion_pass",elt.QuantumDiffusion,"auto"),("energyloss_pass",elt.EnergyLoss,"auto"),("simplequantdiff_pass",elt.SimpleQuantDiff,"auto"),("simpleradiation_pass",elt.SimpleRadiation,"auto"),),True:(("cavity_pass",elt.RFCavity,"auto"),("dipole_pass",elt.Dipole,"auto"),("quadrupole_pass",elt.Quadrupole,"auto"),("wiggler_pass",elt.Wiggler,"auto"),("sextupole_pass",elt.Sextupole,None),("octupole_pass",elt.Octupole,None),("multipole_pass",elt.Multipole,None),("collective_pass",elt.Collective,"auto"),("diffusion_pass",elt.QuantumDiffusion,"auto"),("energyloss_pass",elt.EnergyLoss,"auto"),("simplequantdiff_pass",elt.SimpleQuantDiff,"auto"),("simpleradiation_pass",elt.SimpleRadiation,"auto"),),}# Don't warn on floating-point errorsnp.seterr(divide="ignore",invalid="ignore")warnings.filterwarnings("always",category=AtWarning,module=__name__)# noinspection PyAttributeOutsideInit
[docs]classLattice(list):"""Lattice object. An AT lattice is a sequence of AT elements. In addition to :py:class:`list` methods, :py:class:`Lattice` supports `extended indexing <https://numpy.org/doc/stable/user/basics.indexing.html#advanced-indexing>`_ as a numpy :py:obj:`~numpy.ndarray`. """# Attributes displayed:_disp_attributes=("name","energy","particle","periodicity","harmonic_number","beam_current","nbunch",)# excluded attributes_excluded_attributes=("nbunch",)# Attributes propagated in copies:_std_attributes=("name","_energy","_particle","periodicity","_cell_harmnumber","_radiation","_beam_current","_fillpattern",)# noinspection PyUnusedLocaldef__init__(self,*args,iterator:Filter=None,scan:bool=False,**kwargs):# noinspection PyUnresolvedReferences""" Lattice(elements: Iterable[Element], **params) Lattice(filter, [filter, ...,] iterator=iter, **params) Parameters: elements: iterable of Element objects. iter: function called as :pycode:`iter(params, *args)`. It must return an iterable of :py:class:`.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 (str): Name of the lattice, Default: "". energy (float): Energy of the beam. periodicity (int): Number of periods. Default: 1. If <= 0, it will be deduced from the sum of bending angles. particle (Particle | str): circulating particle. May be "relativistic", "electron", "positron", "proton", "antiproton", "posmuon", "negmuon" or a Particle object. Default: "relativistic". iterator: custom iterator (see below). Default :py:obj:`None`. *: all other keywords will be set as attributes of the Lattice object. beam_current: Total current in the beam, used for collective effects [A]. An iterator ``it`` is called as :pycode:`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 :pycode:`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 :pycode:`Lattice.beam_current = current` The filling pattern and beam current are used by collective effects passmethods. Examples: Chaining iterators (taken from ``load_mat``): >>> ring = Lattice(ringparam_filter, matfile_generator, filename, iterator=params_filter, **params) ``matfile_generator(params, filename)`` opens filename and generates AT elements for each cell of the Matlab cell array representing the lattice, ``ringparam_filter(params, matfile_generator, *args)`` runs through ``matfile_generator(params, *args)``, looks for RingParam elements, fills params with their information and discards them, ``params_filter(params, ringparam_filter, *args)`` runs through ``ringparam_filter(params, *args)``, looks for energy and periodicity if not yet defined. """ifiteratorisNone:(arg1,)=argsor[[]]# accept 0 or 1 argumentifisinstance(arg1,Lattice):elems=lattice_filter(kwargs,arg1)else:elems=params_filter(kwargs,type_filter,arg1)else:elems=iterator(kwargs,*args)super().__init__(elems)# removing excluded attributesforattrinself._excluded_attributes:kwargs.pop(attr,None)# set default valueskwargs.setdefault("name","")periodicity=kwargs.setdefault("periodicity",1)kwargs.setdefault("particle",kwargs.pop("_particle",Particle()))kwargs.setdefault("beam_current",kwargs.pop("_beam_current",0.0))# dummy initialization in case the harmonic number is not therekwargs.setdefault("_fillpattern",np.ones(1))# Remove temporary keywordsfrequency:float|None=kwargs.pop("_frequency",None)cell_length:float|None=kwargs.pop("_length",None)cell_h=kwargs.pop("cell_harmnumber",kwargs.pop("_cell_harmnumber",math.nan))ring_h=kwargs.pop("harmonic_number",cell_h*periodicity)energy=kwargs.setdefault("energy",kwargs.pop("_energy",None))ifenergyisNone:raiseAtError("Lattice energy is not defined")# set attributesself.update(kwargs)# Setting the harmonic number is delayed to have self.beta availableifnot(frequencyisNoneorfrequency==0.0):rev=self.beta*clight/cell_lengthself._cell_harmnumber=int(round(frequency/rev))try:fp=kwargs.pop("_fillpattern",np.ones(1))self.set_fillpattern(bunches=fp)exceptAssertionError:self.set_fillpattern()elifnot((ring_hisNone)ormath.isnan(ring_h)):self._cell_harmnumber=ring_h/periodicitydef__getitem__(self,key):try:# Integerreturnsuper().__getitem__(key.__index__())except(AttributeError,TypeError):ifisinstance(key,slice):# Slicerg=range(*key.indices(len(self)))else:# Array of integers or booleanrg=get_uint32_index(self,key,endpoint=False)returnLattice(elem_generator,(super(Lattice,self).__getitem__(i)foriinrg),iterator=self.attrs_filter,)def__setitem__(self,key,values):try:# Integer or slicesuper().__setitem__(key,values)exceptTypeError:# Array of integers or booleanrg=get_uint32_index(self,key,endpoint=False)fori,vinzip(*np.broadcast_arrays(rg,values)):super().__setitem__(i,v)def__delitem__(self,key):try:# Integer or slicesuper().__delitem__(key)exceptTypeError:# Array of integers or booleanrg=get_uint32_index(self,key,endpoint=False)foriinreversed(rg):super().__delitem__(i)def__repr__(self):ks=", ".join(f"{k}={v!r}"fork,vinself.attrs.items())returnf"Lattice({super().__repr__()}, {ks})"def__str__(self):ks=", ".join(f"{k}={v!r}"fork,vinself.attrs.items())returnf"Lattice(<{len(self)} elements>, {ks})"def__add__(self,elems):"""Add elems, an iterable of AT elements, to the lattice"""returnself.concatenate(elems,copy=True)def__iadd__(self,elems):self.concatenate(elems,copy=False)returnselfdef__mul__(self,n):returnself.repeat(n)def_addition_filter(self,elems:Iterable[Element],copy_elements=False):cavities=[]length=0.0params={}forelemintype_filter(params,elems):ifisinstance(elem,elt.RFCavity):cavities.append(elem)elem.Energy=self._energyelifelem.PassMethod.endswith("RadPass"):elem.Energy=self._energyelifhasattr(elem,"Energy"):delelem.Energyelifhasattr(elem,"_turnhistory"):elem.clear_history(self)length+=getattr(elem,"Length",0.0)ifcopy_elements:yieldelem.deepcopy()else:yieldelemifcavitiesandnothasattr(self,"_cell_harmnumber"):cavities.sort(key=lambdael:el.Frequency)try:self._cell_harmnumber=cavities[0].HarmNumberexceptAttributeError:length+=self.get_s_pos(len(self))[0]rev=self.beta*clight/lengthfrequency=cavities[0].Frequencyself._cell_harmnumber=int(round(frequency/rev))self._radiation|=params.pop("_radiation")
[docs]definsert(self,idx:SupportsIndex,elem:Element,copy_elements=False):r"""This method allow to insert an AT element in the lattice. Parameters: idx (SupportsIndex): index at which the lement is inserted elem (Element): AT element to be inserted in the lattice copy_elements(bool): Default :py:obj:`True`. If :py:obj:`True` a deep copy of elem is used. """# scan the new element to update itelist=list(# noqa: F841self._addition_filter([elem],copy_elements=copy_elements))super().insert(idx,elem)
[docs]defextend(self,elems:Iterable[Element],copy_elements=False):# noinspection PyUnresolvedReferencesr"""This method adds all the elements of `elems` to the end of the lattice. The behavior is the same as for a :py:obj:`list` Equivalent syntaxes: >>> ring.extend(elems) >>> ring += elems Parameters: elems (Iterable[Element]): Sequence of AT elements to be appended to the lattice copy_elements(bool): Default :py:obj:`True`. If :py:obj:`True` deep copies of each element of elems are used """ifhasattr(self,"_energy"):# When unpickling a Lattice, extend is called before the lattice# is initialized. So skip this.elems=self._addition_filter(elems,copy_elements=copy_elements)super().extend(elems)
[docs]defappend(self,elem:Element,copy_elements=False):# noinspection PyUnresolvedReferencesr"""This method overwrites the inherited method :py:meth:`list.append()`, its behavior is changed, it accepts only AT lattice elements :py:obj:`Element` as input argument. Equivalent syntaxes: >>> ring.append(elem) >>> ring += [elem] Parameters: elem (Element): AT element to be appended to the lattice copy_elements(bool): Default :py:obj:`True`. If :py:obj:`True` a deep copy of elem is used """self.extend([elem],copy_elements=copy_elements)
[docs]defrepeat(self,n:int,copy_elements:bool=True):# noinspection SpellCheckingInspection,PyUnresolvedReferences,PyRedeclarationr"""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`. Equivalent syntaxes: >>> newring = ring.repeat(n) >>> newring = ring * n Parameters: n : number of repetitions copy_elements: If :py:obj:`True`, deep copies of the lattice are used for the repetition. Otherwise, the original elements are repeated in the developed lattice. Returns: newring (Lattice): the new repeated lattice """defcopy_fun(elem,copy):ifcopy:returnelem.deepcopy()else:returnelemperiodicity=self.periodicityifn!=0andperiodicity>1:nbp=periodicity/nperiodicity=int(round(nbp))ifabs(periodicity-nbp)>_TWO_PI_ERROR:warnings.warn(AtWarning(f"Non-integer number of cells: {self.periodicity}/{n}."" Periodicity set to 1"),stacklevel=1,)periodicity=1hdict={"periodicity":periodicity}try:hdict.update(harmonic_number=self.cell_harmnumber*n*periodicity)exceptAttributeError:passelems=(copy_fun(el,copy_elements)for_inrange(n)forelinself)returnLattice(elem_generator,elems,iterator=self.attrs_filter,**hdict)
[docs]defconcatenate(self,*lattices:Iterable[Element],copy_elements=False,copy=False):# noinspection PyUnresolvedReferences,SpellCheckingInspection,PyRedeclaration"""Concatenate several `Iterable[Element]` with the lattice Equivalent syntaxes: >>> newring = ring.concatenate(r1, r2, r3, copy=True) >>> newring = ring + r1 + r2 + r3 >>> ring.concatenate(r1, copy=False) >>> ring += r1 Parameters: lattices: :py:obj:`Iterables[Element]` to be concatenanted to the Lattice, several lattices are allowed (see example) copy_elements(bool): Default :py:obj:`False`. If :py:obj:`True` deepcopies of the elements of lattices are used copy(bool): Default :py:obj:`False`. If :py:obj:`True` the lattice is modified in place. Oterwise a new Lattice object is returned Returns: lattice(Lattice): concatenated Lattice, if `copy==True` the new lattice object is returned otherwise None """ifcopy:lattice=Lattice(self)else:lattice=selfforlatinlattices:lattice.extend(lat,copy_elements=copy_elements)returnlatticeifcopyelseNone
[docs]defreverse(self,copy=False):# noinspection PyUnresolvedReferencesr"""Reverse the order of the lattice and swapt the faces of elements. Alignment errors are not swapped Parameters: copy(bool): Default :py:obj:`False`. If :py:obj:`True` the lattice is modified in place. Oterwise a new Lattice object is returned Returns: lattice(Lattice): reversed Lattice, if `copy==True` the new lattice object is returned otherwise None Example: >>> newring = ring.reverse(copy=True) """elems=(el.swap_faces(copy=True)forelinreversed(self))ifcopy:returnLattice(elem_generator,elems,iterator=self.attrs_filter)else:reversed_list=list(elems)self[:]=reversed_list
[docs]defdevelop(self,copy_elements:bool=True)->Lattice:"""Develop a periodical lattice by repeating its elements *self.periodicity* times Parameters: copy_elements: If :py:obj:`True`, deep copies of the elements are used for the repetition. Otherwise, the original elements are repeated in the developed lattice. Returns: newlattice: The developed lattice """returnself.repeat(self.periodicity,copy_elements=copy_elements)
@propertydefattrs(self)->dict:"""Dictionary of lattice attributes"""defextattr(d):forkinself._disp_attributes:d.pop(k,None)yieldk,getattr(self,k,None)vrs=vars(self).copy()# Standard attributesres={k:vfork,vinextattr(vrs)ifvisnotNone}# Custom attributesres.update((k,v)fork,vinvrs.items()ifnotk.startswith("_"))returnres
[docs]defrotate(self,n:int)->Lattice:"""Return a new lattice rotated left by n elements"""iflen(self)==0:returnself.copy()n=n%len(self)# works for n<0returnself[n:]+self[:n]
[docs]defupdate(self,*args,**kwargs)->None:""" update(**kwargs) update(mapping, **kwargs) update(iterable, **kwargs) Update the lattice attributes with the given values """attrs=dict(*args,**kwargs)forkey,valueinattrs.items():setattr(self,key,value)
[docs]defcopy(self)->Lattice:"""Returns a shallow copy of the lattice"""returncopy.copy(self)
[docs]defdeepcopy(self)->Lattice:"""Returns a deep copy of the lattice"""returncopy.deepcopy(self)
[docs]defslice_elements(self,refpts:Refpts,slices:int=1)->Lattice:"""Create a new lattice by slicing the elements at refpts Parameters: refpts: Element selector slices: Number of slices in the specified range. Ignored if size is specified. Default: no slicing Returns: newring: New Lattice object """defslice_generator(_):check=get_bool_index(self,refpts)forel,okinzip(self,check):ifokand(slices>1):frac=np.ones(slices)/slicesyield fromel.divide(frac)else:yieldelreturnLattice(slice_generator,iterator=self.attrs_filter)
[docs]defslice(self,size:float|None=None,slices:int|None=1)->Lattice:"""Create a new lattice by slicing the range of interest into small elements Keyword arguments: size=None: Length of a slice. Default: computed from the range and number of points: ``size = (s_max-s_min)/slices``. slices=1: Number of slices in the specified range. Ignored if size is specified. Default: no slicing Returns: newring: New Lattice object """defslice_generator(_,ibeg,iend):yield fromself[:ibeg]foreleminself[ibeg:iend]:nslices=int(math.ceil(elem.Length/size))ifnslices>1:frac=np.ones(nslices)/nslicesyield fromelem.divide(frac)else:yieldelemyield fromself[iend:]s_range=self.s_rangeifsizeisNone:smin,smax=s_rangesize=(smax-smin)/slicesi_range=self.i_rangereturnLattice(slice_generator,i_range[0],i_range[-1],iterator=self.attrs_filter,s_range=s_range,)
[docs]defattrs_filter(self,params:dict,elem_filter:Filter,*args)->Iterable[Element]:"""Filter function which duplicates the lattice attributes Parameters: params: Dictionary of Lattice attributes elem_filter: Element filter *args: Arguments provided to *elem_filter* """forkeyinself._std_attributes:try:params.setdefault(key,getattr(self,key))exceptAttributeError:passreturnelem_filter(params,*args)
@propertydefs_range(self)->None|tuple[float,float]:"""Range of interest: (s_min, s_max). :py:obj:`None` means the full cell."""try:returnself._s_rangeexceptAttributeError:self.s_range=Nonereturnself._s_range# noinspection PyAttributeOutsideInit@s_range.setterdefs_range(self,value:tuple[float,float]|None):spos=self.get_s_pos(range(len(self)+1))ifvalueisNone:value=(0.0,spos[-1])ok=np.flatnonzero(np.logical_and(spos>value[0],spos<value[1]))iflen(ok)>0:i1=max(ok[0]-1,0)i2=min(ok[-1]+1,len(self))self._i_range=range(i1,i2+1)else:self._i_range=range(0,0)self._s_range=value@propertydefi_range(self)->Uint32Refpts:"""Range of elements inside the range of interest"""try:i_range=self._i_rangeexceptAttributeError:self.s_range=Nonei_range=self._i_rangereturnget_uint32_index(self,i_range)@propertydefenergy(self)->float:"""Lattice energy"""returnself._energy@energy.setterdefenergy(self,energy:float):# Set the Energy attribute of radiating elementsforeleminself:ifisinstance(elem,elt.RFCavity)orelem.PassMethod.endswith("RadPass"):elem.Energy=energy# Set the energy attribute of the Lattice# Use a numpy scalar to allow division by zero# self._energy = np.array(energy, dtype=float)self._energy=energy@propertydefcell_length(self)->float:"""Cell length (1 cell) [m] See Also: :py:meth:`circumference`. """returnself.get_s_pos(len(self))[0]@propertydefcell_revolution_frequency(self)->float:"""Revolution frequency for 1 cell [Hz] See Also: :py:meth:`revolution_frequency`. """beta=self.betareturnbeta*clight/self.cell_length@propertydefcell_harmnumber(self)->float:"""Harmonic number per cell See Also: :py:meth:`harmonic_number`. """try:returnself._cell_harmnumberexceptAttributeErrorasexc:raiseAttributeError("harmonic_number undefined: No cavity found in the lattice")fromexc@propertydefcircumference(self)->float:"""Ring circumference (full ring) [m] See Also: :py:meth:`cell_length`. """returnself.periodicity*self.get_s_pos(len(self))[0]@propertydefrevolution_frequency(self)->float:"""Revolution frequency (full ring) [Hz] See Also: :py:meth:`cell_revolution_frequency`. """beta=self.betareturnbeta*clight/self.circumference@propertydefparticle(self)->Particle:"""Circulating particle"""returnself._particle@particle.setterdefparticle(self,particle:str|Particle):ifisinstance(particle,Particle):self._particle=particleelse:self._particle=Particle(particle)
[docs]defset_wake_turnhistory(self):"""Function to reset the shape of the turn history in collective elements based on the number of slices, turns and bunches """foreinself:ife.is_collective:e.clear_history(ring=self)
[docs]defset_fillpattern(self,bunches:int|np.ndarray=1):"""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). """ifisinstance(bunches,int):ifself.harmonic_number%bunches==0:fp=np.zeros(self.harmonic_number)fp[::int(self.harmonic_number/bunches)]=1else:raiseAtError("Harmonic number has to be a multiple of the scalar input bunches")elifnp.isscalar(bunches):raiseAtError("Scalar input for bunches must be an integer")else:bunches=bunches.astype(dtype=float,casting="safe",copy=False)assert(len(bunches)==self.harmonic_number),f"bunches array input has to be of shape ({self.harmonic_number},)"assertnp.all(bunches>=0.0),"bunches array can contain only positive numbers"fp=bunchesself._fillpattern=fp/np.sum(fp)self.set_wake_turnhistory()
@propertydeffillpattern(self)->np.ndarray:"""Filling pattern describing the bunch relative amplitudes such that ``sum(fillpattern)=1`` """returnself._fillpattern@fillpattern.setterdeffillpattern(self,value):"""Filling pattern describing the bunch relative amplitudes such that ``sum(fillpattern)=1``. Calls the function ``Lattice.set_fillpattern``. """self.set_fillpattern(value)@propertydefbunch_list(self)->np.ndarray:"""Indices of filled bunches"""returnnp.flatnonzero(self._fillpattern)
beam_current=property(get_beam_current,set_beam_current)@propertydefbunch_currents(self)->np.ndarray:"""Bunch currents [A]"""returnself.beam_current*self._fillpattern[self._fillpattern>0]@propertydefbunch_spos(self)->np.ndarray:"""Bunch position around the ring [m]"""try:circ=self.beta*clight*self.harmonic_number/self.rf_frequencyexceptAtError:circ=self.circumferencebs=circ/len(self._fillpattern)allpos=bs*np.arange(len(self._fillpattern))returnallpos[self._fillpattern>0]@propertydefnbunch(self)->int:"""Number of bunches"""# cast to int required from numpy version 2.3.0returnint(np.count_nonzero(self._fillpattern))@propertydefharmonic_number(self)->int:"""Ring harmonic number (full ring) See Also: :py:meth:`cell_harmnumber`. """try:returnint(self.periodicity*self._cell_harmnumber)exceptAttributeErrorasexc:raiseAttributeError("harmonic_number undefined: ""No cavity found in the lattice")fromexc@harmonic_number.setterdefharmonic_number(self,value):cell_h=float(value)/self.periodicity# check on ringifvalue-round(value)!=0:raiseAtError(f"harmonic number ({value}) must be integer")# check on cell# if cell_h-round(cell_h) != 0:# raise AtError('harmonic number ({}) must be a multiple of {}'# .format(value, int(self.periodicity)))self._cell_harmnumber=cell_hiflen(self._fillpattern)!=value:warnings.warn(AtWarning("Harmonic number changed, resetting fillpattern to ""default (single bunch)."),stacklevel=1,)self.set_fillpattern()@propertydefgamma(self)->float:r"""Relativistic :math:`\gamma` of the particles"""rest_energy=self.particle.rest_energyifrest_energy==0.0:rest_energy=e_massreturnfloat(self.energy/rest_energy)@propertydefbeta(self)->float:r"""Relativistic :math:`\beta` of the particles"""rest_energy=self.particle.rest_energyifrest_energy==0.0:return1.0else:gamma=float(self.energy/rest_energy)returnmath.sqrt(1.0-1.0/gamma/gamma)# noinspection PyPep8Naming@propertydefBRho(self)->float:"""Magnetic rigidity [T.m]"""rest_energy=self.particle.rest_energyifrest_energy==0.0:rest_energy=e_massreturnmath.sqrt(self.energy**2-rest_energy**2)/clight@propertydefis_6d(self)->bool:"""True if at least one element modifies the beam momentum See Also: :py:meth:`enable_6d`, :py:meth:`disable_6d`. """try:returnself._radiationexceptAttributeError:radiate=Falseforeleminself:ifelem.longt_motion:radiate=Truebreak# noinspection PyAttributeOutsideInitself._radiation=radiatereturnradiate@propertydefis_collective(self)->bool:""":py:obj:`True` if any element involves collective effects"""foreleminself:ifelem.is_collective:returnTruereturnFalse@propertydefhas_cavity(self)->bool:""":py:obj:`True` if the lattice contains an active :py:class:`RFCavity` """foreleminself:ifelem.PassMethod.endswith("CavityPass"):returnTruereturnFalse# noinspection PyShadowingNames
[docs]defmodify_elements(self,elem_modify:Callable,copy:bool|None=True,**kwargs):"""Modify selected elements, in-place or in a lattice copy Parameters: elem_modify : element selection function. If ``elem_modify(elem)`` returns :py:obj:`None`, the element is unchanged. Otherwise, ``elem_modify(elem)`` must return a dictionary of attribute name and values, to be set to elem. copy: 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 :py:obj:`True`, :py:obj:`None` if copy is :py:obj:`False` Keyword Arguments: """deflattice_modify(**kws):"""Modifies the Lattice with elements modified by elem_modify"""radiate=Falseforeleminself:attrs=elem_modify(elem)ifattrsisnotNone:elem.update(attrs)ifelem.longt_motion:radiate=Trueself._radiation=radiateself.update(kws)deflattice_copy(params):"""Custom iterator for the creation of a new lattice"""radiate=Falseforeleminself:attrs=elem_modify(elem)ifattrsisnotNone:elem=elem.copy()elem.update(attrs)ifelem.longt_motion:radiate=Trueyieldelemparams["_radiation"]=radiateifcopy:returnLattice(lattice_copy,iterator=self.attrs_filter,**kwargs)else:lattice_modify(**kwargs)
def_set_6d(self,enable:bool,*args,**kwargs):"""Set the lattice radiation state"""deflattice_modify():"""Modifies the Lattice in place"""radiate=Falseforeleminself:new_pass=getpass(elem)ifnew_pass:elem.set_longt_motion(enable,new_pass=new_pass,**vargs)ifelem.longt_motion:radiate=Trueself._radiation=radiateself.update(kwargs)deflattice_copy(params):"""Custom iterator for the creation of a new lattice"""radiate=Falseforeleminself:new_pass=getpass(elem)ifnew_pass:elem=elem.set_longt_motion(enable,new_pass=new_pass,copy=True,**vargs)ifelem.longt_motion:radiate=Trueyieldelemparams["_radiation"]=radiatecp=kwargs.pop("copy",False)iflen(args)>0:defgetpass(elem):return"auto"ifisinstance(elem,args)elseNoneifnotall(issubclass(cl,elt.LongtMotion)forclinargs):raiseTypeError("All arguments must be subclasses of 'LongtMotion'")iflen(kwargs)>0:raiseAtError("No keyword is allowed in this mode")else:defgetpass(elem):foreltype,psminpass_table:ifisinstance(elem,eltype):returnpsmreturnNonedefpassm(key,eltype,def_pass):ifallset:def_pass=all_passreturneltype,kwargs.pop(key,def_pass)# Look for global defaultstry:all_pass=kwargs.pop("all_pass")exceptKeyError:allset=Falseelse:allset=True# Build table of PassMethodspass_table=[passm(*defs)fordefsin_DEFAULT_PASS[enable]]vargs={"energy":self.energy}ifenableelse{}ifcp:returnLattice(lattice_copy,iterator=self.attrs_filter,**kwargs)else:lattice_modify()# noinspection PyShadowingNames,PyIncorrectDocstring
[docs]defenable_6d(self,*args,**kwargs)->Lattice|None:# noinspection PyUnresolvedReferencesr""" enable_6d(elem_class[, elem_class]..., copy=False) enable_6d(cavity_pass='auto'[, dipole_pass='auto']..., copy=False) 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: :py: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: * :py:class:`.LongtMotion`: all elements possibly acting on momentum, * :py:class:`.Radiative`: default radiative elements: :py:class:`.Dipole`, :py:class:`.Quadrupole`, :py:class:`.Wiggler`, * :py:class:`.Collective`: all elements dealing with collective effects. The default PassMethod conversion is used, as with the ``'auto'`` keyword value.. **No keyword except** ``copy`` **is allowed in this syntax.** **Syntax using keyword arguments:** Keyword arguments: all_pass: PassMethod overloading the default values for all elements (:py:obj:`None` or 'auto') cavity_pass='auto': PassMethod set on cavities dipole_pass='auto': PassMethod set on dipoles quadrupole_pass='auto': PassMethod set on quadrupoles wiggler_pass='auto': PassMethod set on wigglers sextupole_pass=None: PassMethod set on sextupoles octupole_pass=None: PassMethod set on octupoles multipole_pass=None : PassMethod set on higher order multipoles collective_pass='auto': PassMethod set on collective effect elements (:py:class:`.WakeElement`,...) diffusion_pass='auto': PassMethod set on :py:class:`.QuantumDiffusion` copy=False: If :py:obj:`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. For PassMethod names, the convention is: * :py:obj:`None`: No change * 'auto': Use the default conversion (replace \*Pass by \*RadPass) * anything else: set as the new PassMethod Examples: >>> ring.enable_6d() Modify `ring` in-place, turn cavities ON, turn synchrotron radiation ON in Dipoles and Quadrupoles, turn collective effects ON. >>> ring.enable_6d(at.RFCavity, at.Radiative) Modify `ring` in-place, turn cavities ON, turn synchrotron radiation ON in dipoles and quadupoles. >>> newring = ring.enable_6d(at.Collective, copy=True) Returns a new lattice with collective effects turned ON and nothing else changed >>> newring = ring.enable_6d(all_pass=None, collective_pass="auto", copy=True) Same as the previous example, using the keyword syntax. See Also: :py:meth:`disable_6d`, :py:attr:`is_6d`. """returnself._set_6d(True,*args,**kwargs)
[docs]defdisable_6d(self,*args,**kwargs)->Lattice|None:# noinspection PyUnresolvedReferencesr""" disable_6d(elem_class[, elem_class]... , copy=False) disable_6d(cavity_pass='auto'[, dipole_pass='auto']..., copy=False) Turn longitudinal motion off. By default, remove all longitudinal motion. Modify the lattice in-place or creates a new lattice, depending on the ``copy`` keyword argument. **Syntax using positional arguments:** Parameters: elem_class: :py:class:`.LongtMotion` subclass. Longitudinal motion is turned off for elements which are instances of any given ``elem_class``. In adition to single element classes, a few grouping classes are available: * :py:class:`.LongtMotion`: all elements possibly acting on momentum, * :py:class:`.Radiative`: default radiative elements: :py:class:`.Dipole`, :py:class:`.Quadrupole`, :py:class:`.Wiggler`, * :py:class:`.Collective`: all elements dealing with collective effects. The default PassMethod conversion is used, as with the ``'auto'`` keyword value. **No keyword except** ``copy`` **is allowed in this syntax.** **Syntax using keyword arguments:** Keyword arguments: all_pass: PassMethod overloading the default values for all elements (:py:obj:`None` or 'auto') cavity_pass='auto': PassMethod set on cavities dipole_pass='auto': PassMethod set on dipoles quadrupole_pass=auto: PassMethod set on quadrupoles wiggler_pass='auto': PassMethod set on wigglers sextupole_pass='auto': PassMethod set on sextupoles octupole_pass='auto': PassMethod set on octupoles multipole_pass='auto': PassMethod set on higher order multipoles collective_pass='auto': PassMethod set on collective effect elements (:py:class:`.WakeElement`,...) diffusion_pass='auto': PassMethod set on :py:class:`.QuantumDiffusion` copy=False: If :py:obj:`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. For PassMethod names, the convention is: * :py:obj:`None`: no change * 'auto': Use the default conversion (replace \*RadPass by \*Pass) * anything else: set as the new PassMethod Examples: >>> ring.disable_6d() Modify `ring` in-place, turn OFF everything affecting the longitudinal momentum. >>> ring.disable_6d(at.RFCavity) Turn cavities OFF (nothing else modified). >>> ring.disable_6d(all_pass=None, cavity_pass="auto") Same as the previous example, but using the keyword syntax. >>> newring = ring.disable_6d(cavity_pass=None, copy=True) Return a new Lattice (shallow copy of `ring`) with everything turned OFF except RF cavities. >>> newring = ring.disable_6d( ... all_pass=None, sextupole_pass="DriftPass", copy=True ... ) Return a new Lattice (shallow copy of `ring`) with sextupoles turned into Drifts (turned off) and everything else unchangedd. See Also: :py:meth:`enable_6d`, :py:attr:`is_6d`. """returnself._set_6d(False,*args,**kwargs)
[docs]defsbreak(self,break_s,break_elems=None,**kwargs):r""" Insert elements at selected locations in the lattice. In place modification not available for this function. 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 Example: >>> newring = ring.sbreak(spos, at.Marker('sbreak')) """defsbreak_iterator(_,itmk):"""Iterate over elements and breaks where necessary"""defnext_mk():"""Extract the next element to insert"""try:returnnext(itmk)exceptStopIteration:returnsys.float_info.max,Nones_end=0.0# get the 1st insertionsmk,mk=next_mk()# skip all insertions at negative break_s, if anywhilesmk<s_end:smk,mk=next_mk()foreleminself:s_end+=elem.Length# loop over all insertions within the elementwhilesmk<s_end:frac=(s_end-smk)/elem.Lengthiffrac<1.0:# breakpoint is within the elementel0,elem=elem.divide([1.0-frac,frac])yieldel0yieldmksmk,mk=next_mk()yieldelem# set default insertionifbreak_elemsisNone:break_elems=elt.Marker("sbreak")break_elems=np.reshape(break_elems,-1)# Check element lengthsifnotall(e.Length==0foreinbreak_elems):warnings.warn(AtWarning("Inserting elements with length!=0 may change the lattice"),stacklevel=2,)# broadcast break_s and break_elems to arrays of same size# and create an iterator over the elements to be insertediter_mk=zip(*np.broadcast_arrays(break_s,break_elems))returnLattice(sbreak_iterator,iter_mk,iterator=self.attrs_filter,**kwargs)
[docs]defreduce(self,**kwargs)->Lattice:"""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 Args: keep (Refpts): Keep the selected elements, even with ``IdentityPass`` PassMethod. Default: keep :py:class:`.Monitor` and :py:class:`.RFCavity` elements """defreduce_filter(_,itelem):try:elem=next(itelem).copy()exceptStopIteration:returnfornxtinitelem:try:elem.merge(nxt)exceptTypeError:yieldelemelem=nxt.copy()yieldelemkp=get_bool_index(self,lambdael:el.PassMethod!="IdentityPass")try:keep=get_bool_index(self,kwargs.pop("keep"))exceptKeyError:keep=get_bool_index(self,checktype((elt.Monitor,elt.RFCavity)))returnLattice(reduce_filter,self.select(kp|keep),iterator=self.attrs_filter,**kwargs)
[docs]defreplace(self,refpts:Refpts,**kwargs)->Lattice:"""Return a shallow copy of the lattice replacing the selected elements by a deep copy Parameters: refpts: element selector """check=get_bool_index(self,refpts)elems=(el.deepcopy()ifokelseelforel,okinzip(self,check))returnLattice(elem_generator,elems,iterator=self.attrs_filter,**kwargs)
# Obsolete methods kept for compatibility
[docs]defradiation_on(self,*args,**kwargs)->Lattice|None:"""Obsolete. Turn longitudinal motion on The function name is misleading, since the function deals with longitudinal motion in general. For this reason **the method is obsolete** and replaced by :py:meth:`.enable_6d` See Also: :py:meth:`.enable_6d` """kwargs.update(zip(("cavity_pass","dipole_pass","quadrupole_pass"),args))returnself._set_6d(True,**kwargs)
[docs]defradiation_off(self,*args,**kwargs)->Lattice|None:"""Obsolete. Turn longitudinal motion off The function name is misleading, since the function deals with longitudinal motion in general. For this reason **the method is obsolete** and replaced by :py:meth:`.disable_6d` See Also: :py:meth:`disable_6d` """kwargs.update(zip(("cavity_pass","dipole_pass","quadrupole_pass"),args))returnself._set_6d(False,**kwargs)
@propertydefradiation(self)->bool:"""Obsolete (see :py:attr:`is_6d` instead) :meta private: """try:returnself._radiationexceptAttributeError:radiate=Falseforeleminself:ifelem.longt_motion:radiate=Truebreak# noinspection PyAttributeOutsideInitself._radiation=radiatereturnradiate
[docs]deflattice_filter(params,lattice):"""Copy lattice parameters and run through all lattice elements Parameters: params: Dictionary of Lattice attributes lattice: Input ``Lattice`` Yields: lattice ``Elements`` """returnlattice.attrs_filter(params,elem_generator,lattice)
# noinspection PyUnusedLocal
[docs]defelem_generator(params,elems:Iterable[Element])->Iterable[Element]:"""Run through all elements without any check Parameters: params: Dictionary of Lattice attributes elems: Iterable of lattice ``Elements`` Yields: lattice ``Elements`` """returnelems
no_filter:Filter=elem_generator# provided for backward compatibility
[docs]deftype_filter(params,elems:Iterable[Element])->Generator[Element,None,None]:"""Run through all elements and check element validity. Analyse elements for radiation state Parameters: params: Dictionary of Lattice attributes elems: Iterable of lattice ``Elements`` Yields: lattice ``Elements`` """radiate=Falseforidx,eleminenumerate(elems):ifisinstance(elem,Element):ifelem.longt_motion:radiate=Trueyieldelemelse:warnings.warn(AtWarning(f"item {idx} ({elem}) is not an AT element: ""ignored"),stacklevel=2,)params["_radiation"]=radiate
[docs]defparams_filter(params,elem_filter:Filter,*args)->Generator[Element,None,None]:"""Run through all elements, looking for energy and periodicity. Remove the Energy attribute of non-radiating elements Parameters: params: Dictionary of Lattice attributes elem_filter: Next ``Elements`` filter args: Arguments forwarded to **elem_filter** Yields: lattice ``Elements`` The following keys in ``params`` are set: * ``_length`` * ``periodicity`` * ``energy`` (optional) * ``_frequency`` (optional) * ``harmonic_number`` (optional) energy is taken from: 1) The params dictionary 2) Cavity elements 3) Any other element periodicity is taken from: 1) The params dictionary 2) if periodicity <= 0, from the sum of the bending angles of magnets """el_energies=[]thetas=[]cavities=[]cell_length=0foreleminelem_filter(params,*args):ifisinstance(elem,elt.RFCavity):cavities.append(elem)elifhasattr(elem,"Energy"):el_energies.append(elem.Energy)delelem.Energyifisinstance(elem,elt.Dipole):thetas.append(elem.BendingAngle)cell_length+=getattr(elem,"Length",0.0)yieldelemparams["_length"]=cell_lengthcav_energies=[el.Energyforelincavitiesifhasattr(el,"Energy")]ifcavities:cavities.sort(key=lambdael:el.Frequency)c0=cavities[0]params.setdefault("harmonic_number",getattr(c0,"HarmNumber",math.nan))# params['_harmnumber'] = getattr(c0, 'HarmNumber', math.nan)params["_frequency"]=getattr(c0,"Frequency",None)if"energy"notinparams:energies=cav_energiesorel_energiesifenergies:# Guess energy from the Energy attribute of the elementsenergy=params.setdefault("energy",max(energies))ifmin(energies)<max(energies):warnings.warn(AtWarning("Inconsistent energy values, "f'"energy" set to {energy}'),stacklevel=2,)ifparams.setdefault("periodicity",1)<=0:# Guess periodicity from the bending angles of the latticetry:nbp=2.0*np.pi/sum(thetas)exceptZeroDivisionError:periodicity=1warnings.warn(AtWarning('No bending in the cell, set "Periodicity" to 1'),stacklevel=2,)else:periodicity=int(round(nbp))ifabs(periodicity-nbp)>_TWO_PI_ERROR:warnings.warn(AtWarning(f"Non-integer number of cells: {nbp} -> {periodicity}"),stacklevel=2,)params["periodicity"]=periodicity