lattice#
Lattice building and modifying
Modules
Conversion of AT lattices to/from other formats |
|
Parameter groups |
|
Create lattice elements |
|
Lattice floor plans |
Classes
Particle definition for AT |
Functions
ELSTR=(ELEM) convert AT element tp pyat |
|
Makes the string representation of an AT element |
|
Get the ring properties if existing |
|
Get the ring properties |
|
Add or modify properties of the lattice |
|
adds a multipole component to an existing polynomial, |
|
ataddrandmpole adds a random multipole component to all elements of type |
|
Divide an element into pieces |
|
makes a new AT element structure from another element, |
|
Generate simplified AT structures |
|
Fit chromaticites by scaling 2 sextupole families |
|
Fit linear tunes by scaling 2 quadrupole families |
|
performs a search on MATLAB cell arrays of structures |
|
retrieves the field values AT cell array of elements |
|
Tries to determine the class of an element |
|
extracts the information about element families and |
|
Insert elements at given locations in a line |
|
will load a field error structure into a ring |
|
load a lattice from a file |
|
Get the fundamental mode cavities |
|
MAKERNDFIELDERRS will create a field error data structure |
|
Remove useless elements from an AT structure |
|
Circularly shift the lattice elements |
|
Insert markers at given s positions in a lattice |
|
Set the RF Cavity with the passmethod RFCavityPass. |
|
sets the field values of MATLAB cell array of structures |
|
sets the misalignment vectors |
|
sets the entrance and exit rotation matrices |
|
set new displacement parameters |
|
Creates a “simple ring” |
|
Creates a line by inserting one or more elements into a base element |
|
sets new rotation parameters |
|
Create a JSON file to store an AT lattice |
|
Creates a .m file to store an AT structure |
|
Creates pyAT lattice from a Matlab lattice |
|
places elements from FAMLIST into cell array THERING |
|
combines adjacent elements that have the same specified pass method |
|
combines adjacent elements that use 4-by-5 PassMethods |
|
() Return the list of field names affecting the element entrance |
|
() Return the list of field names affecting the element exit |
|
performs a search on MATLAB cell arrays of structures |
|
looks for string matches in ‘Tag’ field of AT lattice elements |
|
retrieves the field values MATLAB cell array of structures |
|
|
|
inserts one or more elements into a drift element |
|
tests if an input argument is a valid AT element. |
|
tests if an input argument is a valid AT lattice. |
|
removes a lattice element and merges the two adjacent drift spaces |
|
Move an element |
|
Move fields from one structure to another |
|
removes elements of length 0 from the accelerator lattice |
|
sets the field values of MATLAB cell array of structures |
|
sets the misalignment vectors T1, T2 for elements |
|
sets the ‘Tag’ field in AT lattice elements |
|
sets the entrance and exit misalignment matrixes |
|
inserts an element into a drift space |
|
makes a matrix more symplectic |
- class atparticle#
- Particle definition for ATA particle is defined by its name, rest energy and charge
- at2py()#
- ELSTR=(ELEM) convert AT element tp pyatat2py Creates a pyat element
- at2str(elem))#
- Makes the string representation of an AT elementat2str Creates a string such that EVAL(at2str(elem)) recreates anidentical element.INPUTS1. elem - Elem to writeOUTPUTS1. elsstr - String given the AT constructor of the elementSee also
atwritem()
,atwritepy()
- atCheckRingProperties(ring)#
- Get the ring properties if existingproperties=atCheckRingProperties(ring)Extract data from the RingParam element of the lattice, if present,otherwise return an empty structureSee also
atGetRingProperties()
,atSetRingProperties()
- atGetRingProperties(ring, 'param1', 'param2', ...)#
- Get the ring properties[v1,v2,…]=atGetRingProperties(ring,’param1’,’param2’,…)Extract lattice properties. Extract from the RingParam element of thelattice if present, or from the lattice elements. The number of outputscorresponds to the number of property names in input.RING: Ring structureStandard properties (case independent names):‘FamName’ Name of the lattice‘name’ “ “ “ “‘Energy’ Ring energy [eV]‘Periodicity’ Number of periods to build a full ring‘Particle’ particle (Particle object)‘cavpts’ Location of the main cavities‘beta’ Relativistic beta of the particles‘gamma’ Relativistic gamma of the particles‘rf_frequency’ RF frequency (main cavities) [Hz]‘rf_timelag’ RF timelag (main cavities) [m]‘BRho’ Particle rigidity [T.m]‘mcf’ Momentum compaction factor “alpha”‘slip_factor’ Slip factor “eta”‘is_6d’ Longitudinal motion (cavities, radiating elements, …)‘radiation’ “ “‘has_cavity’ Presence of an active RF cavity (implies ‘radiation’)Properties for the full ring (‘periods’ x cells):‘HarmNumber’ Harmonic number (cell_harmnumber * Periodicity)‘harmonic_number’ “ “‘Circumference’ Ring circumference [m] (cell_length * Periodicity)‘rf_voltage’ RF voltage [V] (cell_rf_voltage * Periodicity)‘revolution_frequency’ Revolution frequency [Hz] (cell_revolution_frequency / Periodicity)Properties for one cell:‘cell_harmnumber’ Harmonic number (for 1 cell)‘cell_length’ Cell length [m]‘cell_rf_voltage’ RF voltage per cell [V] (main cavities)‘cell_revolution_frequency’ Revolution frequency [Hz] (for 1 cell)Custom properties may be added with atSetRingProperties. Custom propertynames are case dependent.Example:>> [energy, beta] = atGetRingProperties(ring,’energy’,’beta’);properties=atGetRingProperties(ring)RING: Ring structurePROPERTIES: Structure with fields:‘FamName’ Name of the lattice‘Energy’ Ring energy [eV]‘Periodicity’ Number of periods to build a full ring‘Particle’ particle (Particle object)‘HarmNumber’ Harmonic number (cell_harmnumber * Periodicity)‘cell_harmnumber’ Harmonic number (for 1 cell)‘cavpts’ Location of the main cavitiesproperties=atGetRingProperties(ring,’all’)The PROPERTIES structure contains all the properties described above,both standard and custom.For fast access, some ring properties are stored in a RingParam elementideally located in the 1st position of the lattice. Without such element,the properties are deduced from the lattice contents. This is much slowerand atGetRingProperties displays a warning indicating how to add theRingParam element:>>ring=atSetRingProperties(ring)See also
atSetRingProperties()
- atSetRingProperties(ring [,key,value]...)#
- Add or modify properties of the latticenewring=atSetRingProperties(ring [,key,value]…)Add or modify the attributes of the RingParam element of the lattice,Insert a new RingParam element if necessaryAvailable properties (property names are case-independent):‘FamName’ Name of the lattice‘name’ “ “ “ “‘Energy’ Ring energy [eV]‘Periodicity’ Number of periods to build a full ring‘Particle’ particle (perticle name or Particle object)‘cavpts’ Location of the main cavities‘rf_frequency’ RF frequency (main cavities) [Hz]. Use ‘nominal’to set the nominal frequency‘rf_timelag’ RF timelag (main cavities) [m]Properties for the full ring (‘periods’ x cells):‘HarmNumber’ Harmonic number (cell_harmnumber * Periodicity)‘harmonic_number’ “ “‘rf_voltage’ RF voltage [V] (cell_rf_voltage * Periodicity)Properties for one cell:‘cell_harmnumber’ Harmonic number (for 1 cell)‘cell_rf_voltage’ RF voltage per cell [V] (main cavities)Additional custom fields may be added. They can be retrieved byatGetRingProperties and are saved in files. Custom field names arecase-dependent.For fast access, some ring properties are stored in a RingParam elementideally located in the 1st position of the lattice. If there is no suchelement, atSetRingProperties will add it.See also
atGetRingProperties()
- ataddmpolecomppoly(polynom, index, strength, radius)#
- adds a multipole component to an existing polynomial,scaling it by the reference multipole value and a radius to set the lengthscale[polynomout] = ataddmpolecomppoly(polynom,index,strength,radius)Polynom = given field polynomial that we want to add a component torefindex = index of magnet: 1 for dipole, 2 for quadrupolenewindex: index of Multipole to addstrength: strength of the multipole component at the given radiusradius: reference radius for defining the absolute strengthB^(N)_(n) = radius^(n-N)*b_n/b_N
- ataddmpoleerrors()#
- ataddrandmpole adds a random multipole component to all elements of type‘type’ where type can be ‘dipole’, ‘quadrupole’, or ‘sextupole’[newring] = ATRANDMPOLE(ring,type,newindex,strength,radius)ring = input ringtype = ‘dipole’, ‘quadrupole’ or ‘sextupole’newindex: index of Multipole to addstrength: strength of the multipole component at the given radiusradius: reference radius for defining the absolute strengthif randflag is set to 1, then the errors will be random, GaussiandistributedThe formula for the added errors isB^(N)_(n) = radius^(n-N)*b_n/b_NIt represents the relative field error to the design value at the ref.radiusFor example, to add a random octupole error of 1e-4 at 25 mm, relative to allquadrupoles:newring =ataddrandmpole(ring,’quadrupole’,4,.1e-4,.025,1);
- atdivelem(elem, frac)#
- Divide an element into piecesline=atdivelem(elem,frac)ELEM: Element to be dividedFRAC: Length of each piece, as a fraction of the total lengthLINE: Cell arrayThe sum of frac elements does not need to be 1, however for bendingmagnets, the length will be divided according to FRAC, but the totalbending angle will be divided according to FRAC/SUM(FRAC) so that thetotal bending angle is kept.Example:>> qf=atquadrupole(‘QF’,0.1,0.5);>> line=atdivelem(qf,[0.5;0.5]); % Split a quadrupole in two halvesOptional arguments:‘KeepAxis’, if present, rotations translations are kept at all slices
- atelem(elem, 'field1', 'value1', 'field2', 'value2', ...)#
- makes a new AT element structure from another element,standard AT type, or a template on file. Existing fields are overriddenwith new values. New fields are added and initialized.newelem = atelem(elem,’field1’,’value1’,’field2’, ‘value2’, …)ELEM can be 1) another element structure2) One of the standard AT types‘drift’‘quadrupole’‘sextupole’‘bend’,’rbend’,’sbend’‘marker’‘corrector’…
- atfastring(ring)#
- Generate simplified AT structuresThe given ring structure is modified so that cavities are moved to thebeginning and the rest of the ring is replaced by a linear 6x6 transfermatrix followed by a non-linear element providing tune shifts withamplitude and momentum.[fastring,fastringrad]=atfastring(ring)RING: original AT structure, with no RF and no radiation.FASTRING: Structure containing unchanged cavities moved to thebeginning, a linear 6x6 matrix and a non-linear elementsimulating linear chromaticities and tune shift withamplitudesFASTRINGRAD: Structure containing unchanged cavities moved to thebeginning, a diffusion element, a linear 6x6 transfermatrix and a non-linear element simulating linearchromaticities and tune shift with amplitudes[fastring,fastringrad]=atfastring(ring,refpts)The ring is split at the specified locations, and each section istransformed in the same way as previously described[fastring,fastringrad]=atfastring(ring,’plot’) plots the tune shifts with amplitude
- atfitchrom(ring, newchrom, sextfamily1, sextfamily2)#
- Fit chromaticites by scaling 2 sextupole familiesnewring = atfitchrom(ring,newchrom,sextfamily1,sextfamily2)newring = atfitchrom(ring,dpp,newchrom,sextfamily1,sextfamily2)RING: Cell arrayDPP: Optional momentum deviation (default 0)NEWCHROM: Desired non-normalized chromaticitiesSEXTFAMILY1: 1st sextupole familySEXTFAMILY2: 2nd sextupole familySEXTFAMILY may be:string: Family namelogical array: mask of selected elements in RINGNumeric array: list of selected elements in RINGCell array: All elements selected by each cellnewring = atfitchrom(ring,…,’dpstep’,dpstep)dpstep is the momentum step applied to compute the chromaticity.The default is the DPStep global option, which defaults to 3.0e-6newring = atfitchrom(ring,…,’hstep’,hstep)hstep is the sextupole strength applied to build the jacobian [m^-3].Default: 1.0e-5newring = atfitchrom(ring,…,’dp’,dp)Specify off-momentum if radiation is OFF (default 0)newring = atfitchrom(ring,…,’dct’,dct)Specify path lengthening if radiation is OFF (default 0)See also
atfittune()
- atfittune(ring, newtunes, quadfamily1, quadfamily2)#
- Fit linear tunes by scaling 2 quadrupole familiesnewring = atfittune(ring,newtunes,quadfamily1,quadfamily2)newring = atfittune(ring,dpp,newtunes,quadfamily1,quadfamily2)RING: Cell arrayDPP: Optional momentum deviation (default 0)NEWTUNES: Desired tune valuesQUADFAMILY1: 1st quadrupole familyQUADFAMILY2: 2nd quadrupole familyQUADFAMILY may be:string: Family namelogical array: mask of selected elements in RINGNumeric array: list of selected elements in RINGCell array: All elements selected by each cellnewring = atfittune(ring,…,’useintegerpart’) With this flag, thefunction fits the tunes to the total values of NEWTUNES, includingthe integer part.With this option the function is substantially slower!newring = atfittune(ring,…,’kstep’,kstep)kstep is the quadrupole strength applied to build the jacobian [m^-2].Default: 1.0e-6newring = atfittune(ring,…,’dp’,dp)Specify off-momentum if radiation is OFF (default 0)newring = atfittune(ring,…,’dct’,dct)Specify path lengthening if radiation is OFF (default 0)See also
atfitchrom()
- atgetcells(ring, 'field')#
- performs a search on MATLAB cell arrays of structuresok = atgetcells(ring, ‘field’)returns indexes of elements that have a field named ‘field’ok = atgetcells(ring, ‘field’, value1…)returns indexes of elements whose field ‘field’is equal to VALUE1, VALUE2, … or VALUEN. Where VALUE can either becharacter strings or a number. If its a character string REGULARexpressions can be used.ok = atgetcells(ring, ‘field’, @testfunction, args…)Uses the user-defined TESTFUNCTION to select array elementsTESTFUNCTION must be of the form:OK=TESTFUNTION(ATELEM,FIELDVALUE,ARGS…)OK is a logical array with the same size as RING, refering to matchingelements in RING
- atgetfieldvalues(ring, 'field')#
- retrieves the field values AT cell array of elementsvalues = atgetfieldvalues(ring,’field’) extracts the values ofthe field ‘field’ in all the elements of RINGvalues = atgetfieldvalues(ring,index,’field’) extracts the values ofthe field ‘field’ in the elements of RING selected by INDEXif RING{I}.FIELD is a numeric scalarVALUES is a length(INDEX) x 1 arrayotherwiseVALUES is a length(INDEX) x 1 cell arrayvalues = atgetfieldvalues(…,’default’,default_value) Uses default_valuesif the required field is not existingMore generally atgetfieldvalues(ring,index,subs1,subs2,…) will callGETFIELD(RING{I},subs1,subs2,…) for I in INDEXExamples:v=atgetfieldvalues(ring,1:10,’polynomb’) is a 10x1 cell arraysuch that V{I}=RING{I}.PolynomB for I=1:10v=atgetfieldvalues(ring(1:10),’polynomb’,{1,2}) is a 10x1 arraysuch that V(I)=RING{I},PolynomB(1,2)
- atguessclass(atelem)#
- Tries to determine the class of an elementatclass=atguessclass(atelem) Tries to determine the class of an elementINPUTS1. elem - AT elementNOTES1. atclass=atguessclass(atelem,’useclass’)By default, atguessclass will default “useless” elements (PolynopmB==0)to ‘Drift’ or ‘Marker’, depending on ‘Length’. When specifying‘UseClass’, it it will preserve the ‘Class’ field for those elements.See also
atwritem()
- atindex(lattice)#
- extracts the information about element families andindexes from AT latticeati=atindex(lattice)returns a srtucture with fields named after element familycontaining an array of their AT indexes;ATI.QF = […]ATI.QD = […];…ati=atindex() Uses the global variable THERINGSee also
atgetcells()
- atinsertelems(line, refpts, frac1, elem1[, frac2, elem2...])#
- Insert elements at given locations in a linenewline=atinsertelems(line,refpts,frac1,elem1[,frac2,elem2…])a new line is created by inserting elements at each location specifiedby REFPTS.LINE: Cell array of structuresREFPTS: Insertion points (index list or logical mask)FRAC: Location of the inserted element ELEM within LINE{REFPTS(i)}0<=FRAC<=1if FRAC = 0, ELEM is inserted before LINE{REFPTS(i)} (no splitting)if FRAC = 1, ELEM is inserted after LINE{REFPTS(i)} (no splitting)if FRAC = NaN, LINE{REFPTS(i)} is replaced by ELEM (no check for identical length)if ELEM = [], nothing is inserted, only the splitting takes placeFRAC and ELEM must be scalars or array of the same size as REFPTSSee also
atsplitelem()
,atdivelem()
- atloadfielderrs()#
- will load a field error structure into a ringfield error structure has the fields:elemindex: element indices in ring to impactNval: reference mpole #, e.g. 2 for Quad, 3 for Sextupolenval: multipole index to changeBval: Value of relative normal coefficientAval: Value of relative skew coefficientradius: reference radius used to compute error
- atloadlattice(filepath)#
- load a lattice from a filelattice=atloadlattice(filepath)Create an AT lattice (cell array of AT elements) from FILEPATHlattice=atloadlattice(filepath,’matkey’,variable_name,…)Give the name of a variable containing the lattice in .mat filescontaining several variables. By default, if the .mat file contains asingle variable, it will be loaded. Otherwise atloadlattice will lookin order for ‘ring’, ‘lattice’ and ‘RING’.lattice=atloadlattice(filepath,key,value…)(Key,value) pairs are added to the RingProperties of the lattice oroverload existing ones.Standard keys include: FamName, Energy, Periodicity, HarmNumber, ParticleCustom properties are allowedAllowed file types:.mat Binary Matlab file. If the file contains several variables,the variable name must be specified using the ‘matkey’ keyword..m Matlab function. The function must output a valid AT structure..json JSON filesee also atwritem, atwritejson
- atmaincavities(ring)#
- Get the fundamental mode cavitiesidcav=atmaincavities(ring)Return a refpts pointing to fundamental mode cavities(cavities with the lowest frequency)
- atmakefielderrstruct()#
- MAKERNDFIELDERRS will create a field error data structureclass is an element class, “Quadrupole”, or “Sextupole”
- atreduce(ring)#
- Remove useless elements from an AT structurenewring=atreduce(ring)Remove elements with PassMethod=’IdentityPass’ and merges adjacentsimilar elementsnewring=atreduce(ring,refpts)When merging similar elements, keep REFPTS intact.[newring,kept]=atreduce(…)Returns the index of kept elements so that NEWRING=OLDRING(KEPT)[newring,kept,newrefpts]=atreduce(ring,refpts)Returns in addition the updated list of reference points.
- atrotatelattice(ring, index)#
- Circularly shift the lattice elementsnewring = atrotatelattice(ring,index)Return a new lattice such that RING(INDEX) is the first element
- atsbreak(ring, spos)#
- Insert markers at given s positions in a lattice[newring,refpts]=atsbreak(ring,spos) build a new lattice with breakpointsRING: AT structureSPOS: Poition of the desired breakpointsNEWRING: Modified AT structureREFPTS: Index of breakpoints
- atsetRFCavity(ring, rfv, radflag, harmnumber, deltafreq)#
- Set the RF Cavity with the passmethod RFCavityPass.RFCavityPass allows to change the energy of the beam changing thefrequency of the cavity.newring = atsetRFCavity(ring, rfv, radflag, harmnumber, deltafreq)sets the synchronous phase of the cavity, the voltage, the harmonicnumber and the PassMethod.INPUTS1. ring - Ring structure2. rfv - RF-voltage in volts3. radflag - 0/1: activat/desactivate radiation (atradon/atradoff)4. HarmNumber - Harmonic number5. DeltaFreq - Frequency shift in HzOUTPUTS1. newring - update ring structure with nw RF parametersNOTES1. All the N cavities will have a voltage rfv/N2. radflag says whether or not we want radiation on, which affectssynchronous phase. If radflag is 0, the function calls atradoff,if it is 1, it calls atradon.3. Cavities in the ring must have the Class RFCavity.4. Normally DeltaFreq should be 0, it’s different from 0 when you want tosimulate a change of energy changing the RF frequency. DeltaFreq is inHz.5. Does not work well for misaligned cavityEXAMPLES1. normal use:newring = atsetRFCavity( ring, 6e6, 1, 992, 0 )2. for off-energy simulations:newring = atsetRFCavity( ring, 6e6, 1, 992, 100 )See also
atsetcavity()
,atradon()
,atradoff()
,atgetU0()
- atsetfieldvalues(ring, ...)#
- sets the field values of MATLAB cell array of structuresNote that the calling syntax must be in the form of assignment:ring = atsetfieldvalues(ring,…)MATLAB does not modify variables that only appear on the righthand side as arguments.newring=atsetfieldvalues(ring,’field’,values)In this mode, the function will set values on all the elements of RINGnewring=atsetfieldvalues(ring,index,’field’,values)In this mode, the function will set values on the elements of RINGspecified by INDEX, given as a list of indices or as a logical masknewring=atsetfieldvalues(ring,’field’,valuestruct)In this mode, the function will set values on the elements of RINGwhose family names are given by the field names of VALUESTRUCTnewring=atsetfieldvalues(ring,ringindex,….,valuestruct)As in the previous mode, the function will set values on the elementsof RING whose family names are given by the field names of VALUESTRUCT.But RINGINDEX=atindex(RING) is provided to avoid multiple computations.Field selection———————————————————NEWRING = ATSETFIELD(RING,’field’,VALUES)For each I=1:length(RING), set RING{I}.FIELD=valueNEWRING = ATSETFIELD(RING,’field’,{M,N},VALUES)For each I=1:length(RING), set RING{I}.FIELD(M,N)=valueMore generally,NEWRING = ATSETFIELD(RING,subs1,subs2,…,VALUES)For each I=1:length(RING), SETFIELD(RING{I},subs1,subs2,…,value)The first dimension of VALUES must be either length(INDEX) or 1 (the valuewill be repeated for each element). For a vector to be repeated, encloseit in a cell array.Value format———————————————————Cell array VALUES—————–Mx1 or 1xM cell array : one cell per element1x1 cell array : cell 1 is affected to all selected elementsCharacter array VALUES———————1xN char array (string) : the string as affected to all selected elementsMxN char array : one row per elementNumeric array VALUES——————–1x1 (scalar) : the value is affected to all selected elementsMx1 (column) : one value per elementMxN (matrix) : one row affected per element. If M==1, the single rowis affected to all selected elementsTo assign column vectors to parameters, use a cell array VALUES, with onecolumn per cell
- atsetshift(ring, elemindex, dx, dy)#
- sets the misalignment vectorsring=atsetshift(ring, elemindex, dx, dy) sets the entrance and exit misalignment vectorsof one element or a group of elements in the globally defined lattice THERING.ELEMINDEX contains indexes of elements to be misalignedDX, DY are displacements of the ELEMENTso the coordinate transformation on the particle at entrance isX -> X-DXY -> Y-DYring=atsetshift(ring, elemindex, dx, dy, ‘relativeshift’)adds the shifts DX and DY to the the previous misalignmentsof the elementsatsetshift(elemindex, dx, dy) Uses the global variable THERINGSee also
atsettilt()
- atsettilt(ring, elemindex, psi)#
- sets the entrance and exit rotation matricesof an element or a group of elements in THERINGring=atsettilt(ring,elemindex, psi)ELEMINDEX contains indexes of elements to be rotatedPSI - angle(s) of rotation in RADIANSPOSITIVE PSI corresponds to a CORKSCREW (right)rotation of the ELEMENT looking in the direction of the beam.(or CORKSCREW, aligned with s-axis) rotation of the ELEMENTThe misalgnment matrixes are stored in fields R1 and R2R1 = [ cos(PSI) sin(PSI); -sin(PSI) cos(PSI) ]R2 = R1’ring=atsettilt(ring,elemindex,psi,’relativetilt’)the rotation is added to the previous misalignmentatsettilt(elemindex, psi) Uses the global variable THERINGSee also
atsetshift()
- atshiftelem(oldelem, dx, dz)#
- set new displacement parametersnewelem=atshiftelem(oldelem,dx,dz)DX: horizontal element displacementDY: Vertical element displacementnewelem=atshiftelem(oldelem,dx,dy,’relativeshift’)the shift is added to the previous misalignment of the element
- atsimplering(..., 'name', name, ...)#
- Creates a “simple ring”A “simple ring” consists of:- A RF cavity- a 6x6 linear transfer map- a simple radiation damping element- a detuning and chromaticity element- a simplified quantum diffusion element which provides the equilibriumemittance**ring=atsimplering**(ENERGY,CIRCUMFERENCE,HARMONIC_NUMBER,HTUNE,VTUNE,…VRF,ALPHAC)ENERGY: nominal energy [eV]CIRCUMFERENCE: ring circumference [m]HARMONIC_NUMBER: harmonic numberHTUNE: horizontal tuneVTUNE: vertical tuneVRF: RF voltage [V]ALPHAC: momentum compaction factorRING: generated simple latticering=atsimplering(…,’name’,name,…)specify the name of the generated lattice. Default: ‘’ring=atsimplering(…,’particle’,particle,…)specify the particle. Default: ‘relativistic’ring=atsimplering(…,’qpx’,xsi_h,’qpy’,xsi_v,…)specify the chromaticites. Default: 0ring=atsimplering(…,’a1’,a1,’a2’,a2,’a3’,a3,…)specify the amplitude detuning coefficients. Default: 0ring=atsimplering(…,’betax’,betax,’betay’,betay,…)specify the beta values at origin [m]. Default: 1ring=atsimplering(…,’alphax’,alphax,’alphay’,alpha,…)specify the alpha values at origin. Default: 0ring=atsimplering(…,’dispx’,dispx,’dispxp’,dispxp,…)specify the horizontal dispersion and derivative at origin. Default: 0ring=atsimplering(…,’dispy’,dispy,’dispyp’,dispyp,…)specify the vertical dispersion and derivative at origin. Default: 0ring=atsimplering(…,’emitx’,emitx,’emity’,emity,…)specify the equilibrium emittance. Default: 0, meaning no quantumdiffusionring=atsimplering(…,’taux’,taux,’tauy’,tauy,’tauz’,tauz,…)specify the damping times. Default: 0, meaning no dampingring=atsimplering(…,’espread’,espread,…)specify the energy spread. Default: 0ring=atsimplering(…,’u0’,u0,…)specify the emergy loss per turn [eV]. Default: 0
- atsplitelem(baseelem, frac1, elem1[, frac2, elem2...])#
- Creates a line by inserting one or more elements into a base elementline=atsplitelem(baseelem,frac1,elem1[,frac2,elem2…])Each inserted element is associated with a location given by 0<=FRAC<=1LINE is a cell array containing the sequence of resulting elementsFRACi may be a scalar or a line vector of locations. ELEMi must be asingle element, or a cell array of elements with the same length as FRACi.if FRAC = 0, the element is inserted before BASEELEM (no splitting)if FRAC = 1, the element is inserted after BASEELEM (no splitting)if ELEMi = [], nothing is inserted, only the splitting takes placeatsplitelem will split BASEELEM in elements with negative lengths ifnecessary. Elements with length 0 are not generated. When splittingdipoles, bending angles are distributed as lengths, and face angles areset to zero at splitting points.Examples:>> dr=atdrift(‘DR’,2); % Insert quadrupoles inside a drift>> qf=atquadrupole(‘QF’,0.1,0.5);>> qd=atquadrupole(‘QD’,0.1,-0.5);>> line=atsplitelem(dr,0.2,qf,0.8,qd);>> mk=atmarker(‘MK’);>> line=atsplitelem(qf,0,mk); % Insert a marker before a quadrupole>> line=atsplitelem(qf,0.5,[]); % Split a quadrupole in two halvesSee also
atinsertelems()
,atdivelem()
- attiltelem(oldelem, psi)#
- sets new rotation parametersnewelem = attiltelem(oldelem,psi)PSI - rotation angle in radiansPOSITIVE ROTS corresponds to a CORKSCREW (right)rotation of the ELEMENT looking in the direction of the beam.(or CORKSCREW, aligned with s-axis) rotation of the ELEMENTThe rotation matrixes are stored in fields R1 and R2R1 = [cos(PSI) sin(PSI); -sin(PSI) cos(PSI)]R2 = R1’newelem=attiltelem(oldelem,psi,’relativetilt’)the rotation is added to the previous misalignment
- atwritejson(ring)#
- Create a JSON file to store an AT latticejs=atwritejson(ring)Return the JSON representation of RING as a character arrayatwritejson(ring, filename)Write the JSON representation of RING to the file FILENAMEatwritejson(ring, …, ‘compact’, true)If compact is true, write a compact JSON file (no linefeeds)see also atloadlattice
- atwritem(ring)#
- Creates a .m file to store an AT structureINPUTS1. ring - Lattice structure (.mat file)2. filname - output file where to write the lattice as a line by line file usingelement constructorsEXAMPLES1. atwritem(ring) % Prints the result in the command window2. atwritem(ring,filename) % Prints the result in a fileSee also
at2str()
- atwritepy(ring)#
- Creates pyAT lattice from a Matlab latticepyring=atwritepy(ring)return the python lattice objectpyring=atwritepy(ring,’file’,filename)In addition, store the python lattice object in FILENAME. The ring canbe retrieved in python with the commands:>>> with open(<filename>,’rb’) as fid:… ring=pickle.load(fid)
- buildlat()#
- places elements from FAMLIST into cell array THERINGin the order given by ineger arry ELISTto be use in Accelerator Toolbox lattice definition files
- combinebypassmethod(lattice, method, keepindex, refindex)#
- combines adjacent elements that have the same specified pass method[newlattice, shiftedkeepindex, shiftedref] = combinebypassmethod(lattice,method,keepindex,refindex)
- combinelinear45(lattice, keepindex, refindex)#
- combines adjacent elements that use 4-by-5 PassMethods[newlattice, shiftedkeepindex, shiftedref] = combinelinear45(lattice,keepindex,refindex)
- entrancefields()#
- () Return the list of field names affecting the element entranceOptional arguments:‘KeepAxis’, if present, rotations translations are excluded from listsee also: exitfields atdivelem
- exitfields()#
- () Return the list of field names affecting the element exitOptional arguments:‘KeepAxis’, if present, rotations translations are excluded from listsee also: entrancefields atdivelem
- findcells(cellarray, 'field')#
- performs a search on MATLAB cell arrays of structuresindex = findcells(cellarray, ‘field’)returns indexes of elements that have a field named ‘field’index = findcells(cellarray, ‘field’, value)returns indexes of elements whose field ‘field’is equal to VALUE1, VALUE2, … or VALUEN. Where VALUE can either becharacter strings or a number. If its a character string REGULARexpressions can be used.Example:findcells(thering,’length’,0, 0.2); % will match elements oflengths 0 and 0.2findcells(thering,’famname’,’sfa’,’sda’);
- findtags(cellarray, matchstr)#
- looks for string matches in ‘Tag’ field of AT lattice elementsindex = findtags(cellarray, matchstr)returns indexes of elements that have a field ‘Tag’whose value is a string exactly matching MATCHSTRor a cell array of strings with one element matching MATCHSTRSee also
findcells()
,settags()
- getcellstruct(cellarray, 'field', index, m, n)#
- retrieves the field values MATLAB cell array of structuresvalues = getcellstruct(cellarray,’field’,index,m,n)values = getcellstruct(cellarray,’field’,index,m) can be usedfor one dimensional vectorsvalues = getcellstruct(cellarray,’field’,index) is the same asgetcellstruct(cellarray,’field’,index,1,1) if the field datais a scalarvalues = getcellstruct(cellarray,’field’,index) is a MATLAB cell arrayof strings if specified fields contain strings.
- insertelem0(lattice, driftindex, splitlength, elemdata)#
- - quick and dirty:inserts element(s) of zero length into a drift spacenewlattice = insertelem0(lattice, driftindex, splitlength, elemdata)
- insertindrift(d, qf, 0.5, qd, 2, qf, 3.5)#
- inserts one or more elements into a drift elementand returns a sequence (cell array) of elements ready to to be usedin AT latticeELEMSEQ = INSERTELEM(DRIFT0, ELEM1, POS1, … ELEMN, POSN)EXAMPLE: FODO cell— 1. Declare elementsD = atelem(‘drift’,’Length’,4.5);QF = atelem(‘quad’,’Length’, 1, ‘K’, 1.234);QD = atelem(‘quad’,’Length’, 1, ‘K’, -2.345);— 2. Insert quadrupoles in the drift;fodocell = insertindrift(d, qf, 0.5, qd, 2, qf, 3.5);See also
splitelem()
- isatelem(elem)#
- tests if an input argument is a valid AT element.A valid AT element is a MATLAB structure with requiredfields ‘Length’, ‘PassMethod’, and a set of data fields,specific to the PassMethod used.[test, errorstr] = isatelem(elem)= isatelem(elem, ‘display’)TEST - test result, 1 = valid AT elementERRORSTR - multi-line error messageSee also
passmethod()
,atelem()
- isatlattice(elem[, 'option1', ...])#
- tests if an input argument is a valid AT lattice.A valid AT lattice is a MATLAB cell array of valid AT elements[test, badlist, errorstr] = isatlattice(elem, [‘option1’,…])Allowed otions:‘display’ - print problems to command window;‘stop’ - return after the first problem is foundTEST - test result, 1 = valid AT elementERRORSTR - multi-line error message
- mergedrift()#
- removes a lattice element and merges the two adjacent drift spacesmergedrift (SPLITPOS) removes an element located at SPLITPOS from the global lattice THERINGsurrounded by two DRIFT spaces. The resulting drift has Length L0 = L1 + LSPLIT + L2;Number of elements in THERING is thus reduced by 2See also
splitdrift()
- mvelem(elempos, dist)#
- Move an elementmvelem(elempos, dist) Move an element located at ELEMPOS in THERINGsurrounded by two DRIFT spaces0 < DIST < LD move downstream-LU < DIST < 0 move upstreamwhere LU and LD - lenths ofupstream and downstrem drift drifts BEFORE!!! the moveNumber of elements in THERING and total length remain the sameSee also
splitdrift()
,mergedrift()
- mvfield(dst, src, fieldnames)#
- Move fields from one structure to another[newdst,newsrc]=mvfield(dst,src,fieldnames)Moves the selected fields from SRC to DSTDST: Destination structureSRC: Source structureFIELDNAMES: Field names to be moved (Default: all fields from SRC)NEWDST: DST structure with fields addedNEWSRC: SRC structure with fields removed
- rmelem0(lattice, elemindex)#
- removes elements of length 0 from the accelerator latticenewlattice = rmelem0(lattice,elemindex)[newlattice, shiftedindex] = rmelem0(lattice,elemindex)The number of elements in the modified lattice isreduced by length(ELEMINDEX)SHIFTEDINDEX points to elements in the NEWLATTICE thatimmediateley followed the removed elements in the original LATTICESee also
splitdrift()
,mergedrift()
- setcellstruct(cellarray, ...)#
- sets the field values of MATLAB cell array of structuresNote that the calling syntax must be in the form of assignment:cellarray = setcellstruct(cellarray,…)MATLAB does not modify variables that only appear on the righthand side as arguments.Numeric data———————————————————cellarray = setcellstruct(cellarray,’field’,index,value,m,n)Sets (M,N) element equal to VALUE when the field data isa matrix. The assigned VALUE may be1. Scalar numeric value to be written to all CELLARRAY elements2. Numeric array of the same length as INDEX arraycellarray = setcellstruct(cellarray,’field’,index,value,m) can be usedfor one dimensional vectorscellarray = setcellstruct(cellarray,’field’,index,value) is the same assetcellstruct(cellarray,’field’,index,1,1) if the field datais a scalarCharacter array——————————————————————–CELLARRAY setcellstruct(cellarray,’field’,index,value) is a MATLABcell array of strings when specified fields contain strings.The assignment VALUE may be1. Character string,2. Character array with the number of rows matching the number ofelements in INDEX array,3. Cell array of strings, with either one element or with the samelength as index.
- setshift(elemindex, dx, dy)#
- sets the misalignment vectors T1, T2 for elementssetshift(elemindex, dx, dy) sets the entrance and exit misalignment vectorsof one element or a group of elements in the globally defined lattice THERING.DX, DY are displacements of the ELEMENTso the coordinate transformation on the particle at entrance isX -> X-DXY -> Y-DYThe elements to be modified are given by ELEMINDEXPrevious stored values are overwritten.See also
setshift()
- settags(lattice, index, tag)#
- sets the ‘Tag’ field in AT lattice elementslattice = settags(lattice, index, tag)INDEX can be integer AT index or a string famly nameTAG is a string tag or a cell array of stringslattice = settags(lattice, index, tag, ‘append’)appends to existing tags
- settilt(elemindex, psi)#
- sets the entrance and exit misalignment matrixesof an element or a group of elements in THERINGPreviously stored values are overwritten.settilt(elemindex, psi)ELEMINDEX contains indexes of elements to be rotatedPSI - angle(s) of rotation in RADIANSPOSITIVE PSI corresponds to a CORKSCREW (right)rotation of the ELEMENT.(or CORKSCREW, aligned with s-axis) rotation of the ELEMENTThe misalgnment matrixes are stored in fields R1 and R2R1 = [ cos(PSI) sin(PSI); -sin(PSI) cos(PSI) ]R2 = R1’NOTES1. This function is deprecated. Use atsettilt insteadSee also
setshift()
,mksrollmat()
- splitdrift(driftpos, split) inserts a marker (zero-length)#
- inserts an element into a drift spacesplitdrift(driftpos, split) inserts a marker (zero-length) elementat distance SPLIT ( 0 < SPLIT < 1) into a drift spacelocated at DRIFTPOS in THERINGsplitdrift(driftpos, split, elemstruccture) inserts a marker (zero-length) elementat distance SPLIT ( 0 < SPLIT < 1) into a drift spacelocated at DRIFTPOS in THERINGNumber of elements in the RING is thus increased by 2SPLIT (controls the position of the splitL1 = L0*SPLITL2 = L0(1-SPLIT)where L0 is the length of the original DRIFTSee also
mergedrift()
- symplectify()#
- makes a matrix more symplecticfollow Healy algorithm as described by McKayBNL-75461-2006-CP