Class CalculationsOnMachines

java.lang.Object
xal.tools.beam.calc.CalculationEngine
xal.tools.beam.calc.CalculationsOnMachines
All Implemented Interfaces:
ISimulationResults, ISimulationResults.ISimEnvResults<TransferMapState>, ISimulationResults.ISimLocResults<TransferMapState>
Direct Known Subclasses:
CalculationsOnRings

Class for performing the calculations expressed in the ISimEnvResults interface in the context of a particle beam system without regard to the particle. That is, only properties of the machine are computed, no attributes or properties of the beam are required for the computations.
Since:
Nov 7, 2013
Author:
Christopher K. Allen
  • Constructor Details

    • CalculationsOnMachines

      public CalculationsOnMachines(Trajectory<TransferMapState> datSim) throws IllegalArgumentException

      Constructor for CalculationsOnMachines. Accepts the TransferMapTrajectory object and extracts the final state and full trajectory transfer map. Quantities that are required for subsequent machine property calculations are also computed, such as phase advance through the trajectory (modulo 2π), entrance position "fixed orbit" (that is, the orbit that is invariant for repeated applications of the linear part of the transfer map (separating the projective part), and the matched envelope when treating the full transfer map as the map for a periodic cell.

      Parameters:
      datSim - the simulation data for the ring, a "transfer map trajectory" object
      Throws:
      IllegalArgumentException - the trajectory does not contain TransferMapState objects
      Since:
      Nov 7, 2013
  • Method Details

    • computeTransferMatrix

      public static PhaseMatrix computeTransferMatrix(TransferMapState state1, TransferMapState state2)
      Convenience method for computing the transfer matrix between two state locations, say S1 and S2. Let s0 be the axis location of the beamline entrance, s1 the location of state S1, and s2 the location of state S2. Each state object Sn contains the transfer matrix Φ(sn,s0) which takes phases coordinates at the beamline entrance to the position of state Sn. The transfer matrix Φ(s2,s1) taking phase coordinates z1 (and covariance matrix σ1) from position s1 to position s2 is then given by

          Φ(s2,s1) = Φ(s2,s0) Φ(s1,s0)-1 ,

      where Φ(s2,s0) is the transfer matrix between the beamline entrance s0 and the position s2 of state S2, and Φ(s1,s0) is the transfer matrix between the beamline entrance s0 and the position s1 of state S1.
      Parameters:
      state1 - trajectory state S1 of starting location s1
      state2 - trajectory state S2 of final location s2
      Returns:
      transfer matrix Φ(s2,s1) between locations s1 and s2
      Since:
      Jun 23, 2014
    • computeTransferMap

      public static PhaseMap computeTransferMap(TransferMapState state1, TransferMapState state2)
      Convenience method for computing the transfer map between two state locations, say S1 and S2. Let s0 be the axis location of the beamline entrance, s1 the location of state S1, and s2 the location of state S2. Each state object Sn contains the transfer map T(sn,s0) which takes phases coordinates at the beamline entrance to the position of state Sn. The transfer map T(s2,s1) taking phase coordinates z1 (and covariance matrix σ1) from position s1 to position s2 is then given by

          T(s2,s1) = T(s2,s0) ∘ T(s1,s0)-1 ,

      where T(s2,s0) is the transfer map between the beamline entrance s0 and the position s2 of state S2, and T(s1,s0) is the transfer map between the beamline entrance s0 and the position s1 of state S1.
      Parameters:
      state1 - trajectory state S1 of starting location s1
      state2 - trajectory state S2 of final location s2
      Returns:
      transfer map T(s2,s1) between locations s1 and s2
      Since:
      Nov 4, 2014
    • getTrajectory

      public Trajectory<TransferMapState> getTrajectory()
      Returns the simulation trajectory from which all the machine properties are computed.
      Returns:
      simulation trajectory from which this object was initialized
      Since:
      Nov 7, 2013
    • getFullTransferMap

      public PhaseMap getFullTransferMap()
      Returns the transfer map of the full machine lattice represented by the associated simulation trajectory.
      Returns:
      the transfer map of the last state of the associated trajectory
      Since:
      Nov 7, 2013
    • getMatchedTwiss

      public Twiss[] getMatchedTwiss()

      Returns the collection of matched Courant-Snyder parameters for each phase plane. The Courant-Snyder parameters describe the matched beam envelopes when treating the trajectory of transfer maps as that for a periodic lattice. There the beam envelopes would have the same characteristics at the beginning and end of the lattice.

      These are computed parameters calculated once at construction time.

      Returns:
      the Courant-Snyder parameters of a beam needed to match it into the transfer maps of the associated trajectory when treated as a periodic lattice
      Since:
      Nov 7, 2013
    • computeTransferMatrix

      public PhaseMatrix computeTransferMatrix(String strIdElemFrom, String strIdElemTo)

      Returns the state response matrix calculated from the front face of elemFrom to the back face of elemTo. This is a convenience wrapper to the real method in the trajectory class

      This method was moved here from EnvelopeTrajectory/EnvelopeProbe where it was eliminated since Trajectory was genericized.

      Parameters:
      strIdElemFrom - String identifying starting lattice element
      strIdElemTo - String identifying ending lattice element
      Returns:
      response matrix from elemFrom to elemTo
      See Also:
      • EnvelopeTrajectory#computeTransferMatrix(String, String)
    • computeCoordinatePosition

      public PhaseVector computeCoordinatePosition(TransferMapState state)

      We return the projective portion of the full-turn transfer map φn : P6P6 where n is the index of the given state Sn. This is the image Δz of the value 0P6R6 × {1}. That is the value Δz = φn(0).

      Recall that the transfer map φ is a PhaseMap object containing a first-order component which is a linear operator on projective space P6. As such, this PhaseMatrix object Φ is embedded in R7×7. The 7th column Δz of Φ is the column of translation operations on a phase vector zP6R6 × {1} since zP6 is represented

          z = (x, x', y, y', z, z', 1)T .

      Thus, the action of Φ can be (loosely) decomposed as

          Φ ⋅ z = Mz + Δz ,

      where MSp(6), the symplectic group. This method returns the component Δz of the transfer map.

      Specified by:
      computeCoordinatePosition in interface ISimulationResults.ISimLocResults<TransferMapState>
      Parameters:
      state - state containing location of returned translation vector
      Returns:
      the translation vector Δz of the given transfer map
      Since:
      Oct 22, 2013
    • computeFixedOrbit

      public PhaseVector computeFixedOrbit(TransferMapState state)

      Get the fixed point at this state location about which betatron oscillations occur.

      Calculate the fixed point solution vector representing the closed orbit at the location of this element. We first attempt to find the fixed point for the full six phase space coordinates. Let Φ denote the one-turn map for a ring. The fixed point zR6×{1} in homogeneous phase space coordinates is that which is invariant under Φ, that is,

          Φz = z .

      This method returns that vector z.

      Recall that the homogeneous transfer matrix Φ for the ring has final row that represents the translation Δ of the particle for the circuit around the ring. The 6×6 sub-matrix of Φ represents the (linear) action of the bending magnetics and quadrupoles and corresponds to the matrix TR6×6 (here T is linear). Thus, we can write the linear operator Φ as the augmented system

           Φ = |T Δ |,   z ≡ |p| ,
               |0 1 |        |1|
       
      where p is the projection of z onto the ambient phase space R6 (without homogeneous the homogeneous coordinate).

      Putting this together we get

          Φz = Tp + Δ = p ,

      to which the solution is

          p = -(T - I)-1Δ

      assuming it exists. The question of solution existence falls upon the resolvent R ≡ (T - I)-1 of T. By inspection we can see that p is defined so long as the eigenvalues of T are located away from 1. In this case the returned value is the augmented vector (p 1)TR6 × {1}.

      When the set of eigenvectors does contain 1, we attempt to find the solution for the transverse phase space. That is, we take vector pR4 and TR4×4 where T = proj4×4 Φ. The returned value is then z = (p 0 0 1)T.

      Specified by:
      computeFixedOrbit in interface ISimulationResults.ISimLocResults<TransferMapState>
      Parameters:
      state - state containing transfer map and location used in these calculations
      Returns:
      fixed point solution (x,x',y,y',z,z',1) for the given phase matrix
      Since:
      Oct 25, 2013
    • computeChromAberration

      public PhaseVector computeChromAberration(TransferMapState state)
      Computes the chromatic aberration for one pass around the ring starting at the given state location, or from the entrance to state position for a linear machine. The returned vector is the displacement from the closed orbit caused by a unit momentum offset (δp = 1). See the documentation in ISimLocResults#computeChromAberration(ProbeState) for a more detailed exposition.
      Specified by:
      computeChromAberration in interface ISimulationResults.ISimLocResults<TransferMapState>
      Parameters:
      state - simulation state where parameters are computed
      Returns:
      the vector Δ of dispersion coefficients
      Since:
      Nov 15, 2013
      See Also:
      • xal.tools.beam.calc.ISimLocResults#computeChromAberration(xal.model.probe.traj.ProbeState)
    • computeTwissParameters

      public Twiss[] computeTwissParameters(TransferMapState state)

      Computes the closed orbit Twiss parameters at this state location representing the matched beam envelope when treating the trajectory as a lattice of periodic cells.

      Returns the array of twiss objects for this state for all three planes.

      Calculates the matched Courant-Snyder parameters for the given period cell transfer matrix and phase advances. When the given transfer matrix is the full-turn matrix for a ring the computed Twiss parameters are the matched envelopes for the ring at that point.

      Let Φ denote the transfer matrix from the machine beginning to the given state location (it is contained in the given state). It is assumed to be the transfer matrix through at least one cell in a periodic lattice. Internally, the array of phase advances {σx, σy, σx} are assumed to be the particle phase advances through the cell for the matched solution. These are computed with the method CalculationEngine.calculatePhaseAdvPerCell(PhaseMatrix).

      The returned Courant-Snyder parameters (α, β, ε) are invariant under the action of the given phase matrix, that is, they are matched. All that is require are α and β since ε specifies the size of the beam envelope. Consequently the returned ε is NaN.

      The following are the calculations used for the Courant-Snyder parameters of a single phase plane:

          α ≡ -ww' = (φ11 - φ22)/(2 sin σ) ,

          β ≡ w2 = φ12/sin σ

      where φij are the elements of the 2×2 diagonal blocks of Φ corresponding the the particular phase plane, the function w is taken from Reiser, and σ is the phase advance through the cell for the particular phase plance.

      Specified by:
      computeTwissParameters in interface ISimulationResults.ISimEnvResults<TransferMapState>
      Parameters:
      state - state containing transfer map and location used in these calculations
      Returns:
      array (twiss-H, twiss-V, twiss-L)
      Since:
      Aug 14, 2013
    • computeBetatronPhase

      public R3 computeBetatronPhase(TransferMapState state)

      This is the phase advance for the given state location.

      Compute and return the particle phase advance from the trajectory beginning to the given state location.

      Internally the method calculates the phase advances using the initial and final Courant-Snyder α and β values for the matched beam at the trajectory beginning and the location of the given state, respectively. These Courant-Snyder parameters are computed as the matched beam envelope at the trajectory beginning and the given state location.

      Let Φ represent the transfer matrix from the initial ring location to the given state location. The computed quantity is the general phase advance of the particle through the transfer matrix Φ, and no special requirements are placed upon Φ (e.g., periodicity). One phase advance is provided for each phase plane, i.e., (σx, σz, σz).

      The definition of phase advance σ is given by

          σ(s) ≡ ∫s [1/β(t)]dt ,

      where β(s) is the Courant-Snyder, envelope function, and the integral is taken along the interval between the initial and final Courant-Snyder parameters.

      The basic relation used to compute σ is the following:

          σ = sin-1 φ12/(β1β2)½ ,

      where φ12 is the element of Φ in the upper right corner of each 2×2 diagonal block, β1 is the initial beta function value (provided) and β2 is the final beta function value (provided).

      Specified by:
      computeBetatronPhase in interface ISimulationResults.ISimEnvResults<TransferMapState>
      Parameters:
      state - state containing transfer map and location used in these calculations
      Returns:
      vector (ψx, ψy, ψx) of phases in radians
      Since:
      Aug 14, 2013
      See Also:
    • computeChromDispersion

      public PhaseVector computeChromDispersion(TransferMapState state)

      Calculates the fixed point (closed orbit) in transverse phase space at the given state location in the presence of dispersion.

      Let the full-turn map a the state location be denoted Φ. The transverse plane dispersion vector Δ is defined

          Δt ≡ -(1/γ2)[dx/dz', dx'/dz', dy/dz', dy'/dz']T .

      It can be identified as the first 4 entries of the 6th column in the transfer matrix Φ. The above vector quantifies the change in the transverse particle phase coordinate position versus the change in particle momentum. The factor -(1/γ2) is needed to convert from longitudinal divergence angle z' used by XAL to momentum δp ≡ Δp/p used in the dispersion definition. Specifically,

          δp ≡ Δp/p = γ2z'

      As such, the above vector can be better described

          Δt ≡ [Δxp, Δx'p, Δyp, Δy'p]T

      explicitly describing the change in transverse phase coordinate for fractional change in momentum δp.

      Since we are only concerned with transverse phase space coordinates, we restrict ourselves to the 4×4 upper diagonal block of Φ, which we denote take T. That is, T = π ⋅ Φ where π : R6×6R4×4 is the projection operator.

      This method finds that point zt ≡ (xt, x't, yt, y't) in transvse phase space that is invariant under the action of the ring for a given momentum spread δp. That is, the particle ends up in the same location each revolution. With a finite momentum spread of δp > 0 we require this require that

          Tzt + δpΔt = zt ,

      which can be written

        zt = δp(T - I)-1Δt ,

      where I is the identity matrix. Dividing both sides by δp yields the final result

        z0ztp = (T - I)-1Δt ,

      which is the returned value of this method. It is normalized by δp so that we can compute the closed orbit for any given momentum spread.

      Specified by:
      computeChromDispersion in interface ISimulationResults.ISimEnvResults<TransferMapState>
      Parameters:
      state - we are calculating the dispersion at this state location
      Returns:
      The closed orbit fixed point z0 for finite dispersion, normalized by momentum spread. Returned as an array [x0,x'0,y0,y'0]/δp
      Since:
      Nov 8, 2013
      See Also:
      • xal.tools.beam.calc.ISimEnvResults#computeChromDispersion(xal.model.probe.traj.ProbeState)
    • calculateFullLatticeMatrixAt

      protected PhaseMatrix calculateFullLatticeMatrixAt(TransferMapState state)

      Calculates and returns the full lattice matrix for the machine at the given state location. Let Sn be the given state object at location sn, and let Tn be the transfer matrix between locations s0 and sn , where s0 is the location of the full transfer matrix Φ0 for this machine (end to end). Then the full turn matrix Φn for the machine at location sn is given by

          Φn = TnΦ0Tn-1 .

      That is, we conjugate the full transfer map for this machine by the transfer map for the given state.

      Parameters:
      state - state object Sn for location sn containing transfer matrix Tn
      Returns:
      the full-turn map Φn at the location sn of the given state
      Since:
      Oct 28, 2013