Class BeamEllipsoid

java.lang.Object
xal.tools.beam.em.BeamEllipsoid

public class BeamEllipsoid extends Object

Encapsulates the properties of a ellipsoidally symmetric distribution of charge in 6D phase space, in particular, the electromagnetic properties. Provides convenience methods for creating arbitrarily oriented ellipsoids and determining their fields, and their effects on beam particles, specifically, in the form of a linear transfer matrix generator.

The ellipsoid is assumed to be moving in the axial (z-axis) direction in the laboratory frame. All relativistic calculations are made with this assumption. In particular, the gamma (relativistic factor) is taken with respect to this axis.

The reference ellipsoid in three-space is represented by the equation

    r'σ-1r = 1

where σ is the 3×3 matrix of spatial moments, r ≡ (x y z) is the position vector in R3, and the prime indicates transposition.

The matrix σ is taken from the 7×7 covariance matrix τ in homogeneous coordinates; it is the matrix formed from the six, second-order spatial moments <xx>, <xy>, <xz>, <yy>, <yz>, <zz> and is arranged as follows:

    | <xx> <xy> <xz> |
σ = | <xy> <yy> <yz> |
    | <xz> <yz> <zz> |


Note that σ must be symmetric and positive definite. Thus, it is diagonalizable with the decomposition

    σ = RΛR'

where the prime indicates transposition, RSO(3) is the orthogonal rotation matrix and Λ is the diagonal matrix of real eigenvalues. Note from the above that these eigenvalues are the squares of the ellipsoid semi-axes.

Author:
Christopher K. Allen
  • Field Summary

    Fields
    Modifier and Type
    Field
    Description
    static final double
    distribution dependence factor - use that for the uniform beam
  • Constructor Summary

    Constructors
    Constructor
    Description
    BeamEllipsoid(double dblGamma, CovarianceMatrix matSigLab)
    Construct a beam charge density ellipsoid described by the phase space correlation matrix matSigLab and relativistic factor gamma.
    BeamEllipsoid(double dblGamma, R3 vec1stMmts, R3 vec2ndMmts)
    Construct a beam charge density ellipsoid where the second spatial moments and axial displacement are given directly.
  • Method Summary

    Modifier and Type
    Method
    Description
    static double[]
    compDefocusConstants(double dblGamma, double[] arrMoments)
    Computes the normalized defocusing constants (squared) for a space charge kick given the second moments of the beam ellipsoid.
    static double[]
    compDefocusConstantsAlaTrace3D(double dblGamma, double[] arrMoments)
    This method is provided as a comparison utility for validation against simulation with Trace3D.
    computeDCScheffMatrix(double dblLen, double dblPerveance)
     
    computeScheffGenerator(double dblPerveance)
    Calculates the transfer matrix generator for space charge effects from this BeamEllipsoid object given the generalized beam three-dimensional perveance dblPerveance.
    computeScheffMatrix(double dblLen, double dblPerveance)
    Compute and return the transfer matrix for space charge effects due to this beam ellipsoid for the given incremental path length dblLen and generalized beam dblPerveance.
    double[]
    Return all the ellipsoid second spatial moments in the stationary beam frame and aligned to the coordinate axes.
    double
    Return the value of the first ellipsoid second spatial moment in the stationary beam frame and aligned to the coordinate axes.
    double
    Return the value of the second ellipsoid second spatial moment in the stationary beam frame and aligned to the coordinate axes.
    double
    Return the value of the third ellipsoid second spatial moment in the stationary beam frame and aligned to the coordinate axes.
    Return the complete transformation from the beam inertial coordinates to the ellipsoid inertial coordinates.
    Return the correlation matrix for the beam in the stationary beam coordinates.
    Return the original correlation matrix for the beam in the laboratory coordinates.
    double[]
    Return all the normalized space-charge defocusing constants for the beam ellipsoid.
    double
    Return the 1st normalized space charge defocusing constant.
    double
    Return the 2nd normalized space charge defocusing constant.
    double
    Return the 3rd normalized space charge defocusing constant.
    double
    Return the relativistic parameter for the ellipsoidal charge distribution.
    Return the complete transformation from the laboratory inertial coordinates to the ellipsoid inertial coordinates.
    Get the Lorentz transform matrix which takes the laboratory coordinates to the beam coordinates.
    Get orthogonal rotation matrix R in SO(7) that rotates the ellipsoid spatial semi-axes onto the spatial coordinate axes.
    double[]
    Return all the ellipsoid semi-axes lengths as an array.
    double
    Return the value of the first ellipsoid semi-axis in the stationary beam frame.
    double
    Return the value of the second ellipsoid semi-axis in the stationary beam frame.
    double
    Return the value of the third ellipsoid semi-axis in the stationary beam frame.
    Return the translation matrix which transforms coordinates in the beam frame to those with the ellipsoid centroid as the coordinate origin.

    Methods inherited from class java.lang.Object

    clone, equals, finalize, getClass, hashCode, notify, notifyAll, toString, wait, wait, wait
  • Field Details

    • CONST_UNIFORM_BEAM

      public static final double CONST_UNIFORM_BEAM
      distribution dependence factor - use that for the uniform beam
  • Constructor Details

    • BeamEllipsoid

      public BeamEllipsoid(double dblGamma, CovarianceMatrix matSigLab)

      Construct a beam charge density ellipsoid described by the phase space correlation matrix matSigLab and relativistic factor gamma.

      The correlation matrix matSigLab is taken in laboratory coordinates in the laboratory frame. The beam is assumed to be moving in the axial (z-axis) direction in the laboratory frame with relativistic factor dblGamma.

      Note that the phase space correlation matrix in homogeneous coordinates contains all moments up to and including second order; this includes the first-order moments that describe the displacement of the ellipsoid from the coordinate origin. Thus, from the 7x7 phase space correlation matrix we extract the displacement vector

          (<x> <y> <z>)

      and the 3×3 covariance matrix σ

          | <x*x> <x*y> <x*z> |   | <x>*<x> <x>*<y> <x>*<z> |
      σ | <y*x> <y*y> <y*z> | - | <y>*<x> <y>*<y> <y>*<z> |
          | <z*x> <z*y> <z*z> |   | <z>*<x> <z>*<y> <z>*<z> |


      to construct the ellipsoidal charge object according to the class documentation

      Parameters:
      dblGamma - relativistic factor
      matSigLab - envelope correlation matrix in homogeneous phase space coordinates
    • BeamEllipsoid

      public BeamEllipsoid(double dblGamma, R3 vec1stMmts, R3 vec2ndMmts)

      Construct a beam charge density ellipsoid where the second spatial moments and axial displacement are given directly. These values are taken in the laboratory coordinates and the beam is assumed to be moving in the axial (z-axis) direction in the laboratory frame with relativistic factor dblGamma.

      The arguments provided include the first-order moments vec1stMmts describing the displacement of the beam ellipsoid from the coordinate origin. This is the vector

          (<x> <y> <z>)

      The argument vec2ndMmts should be the central spatial moments <x2> - <x>2, <y2> - <y>2, <z2> - <z>2, taken from the central (spatial) covariance matrix

          | <x*x> <x*y> <x*z> |   | <x>*<x> <x>*<y> <x>*<z> |
      σ | <y*x> <y*y> <y*z> | - | <y>*<x> <y>*<y> <y>*<z> |
          | <z*x> <z*y> <z*z> |   | <z>*<x> <z>*<y> <z>*<z> |


      In this constructor we assume that the beam ellipsoid is aligned with the laboratory coordinate axes. Thus, no rotations are necessary and the process of computing the space charge effects is expedited.

      NOTES:

      This constructor has not yet been tested and debugged!

      Parameters:
      dblGamma - the relativistic factor
      vec1stMmts - vector (<x> <y> <z>) of first spatial moments
      vec2ndMmts - vector (<x2> - <x>2, <y2> - <y>2, <z2> - <z>2) of second spatial moments
      Since:
      Aug 25, 2011
  • Method Details

    • compDefocusConstants

      public static double[] compDefocusConstants(double dblGamma, double[] arrMoments)

      Computes the normalized defocusing constants (squared) for a space charge kick given the second moments of the beam ellipsoid. The true defocusing constants are obtained by multiplying the result by the factor K*ds, i.e.,

          k2 = kn2K ds

      where k2 is the true defocusing constant, kn2 is the normalized defocusing constant(returned by this method), K is the generalized beam perveance, and ds is the incremental path length over which the kick is being applied.

      The defocusing constants (knx,kny,knz) are given for the Cartesian coordinate system (x,y,z) which is aligned to the ellipsoid semi-axes. If the ellipsoid is rotated with respect to the coordinate system in which the ellipsoid was defined (constructed), then any transfer matrix built from these focusing constants must be rotated into the original coordinate system.

      In the natural coordinates of the ellipsoid the normalized focusing constants kn are related to the space charge defocusing lengths f by the formula

          f = 1/(ds K kn2)

      From an optics standpoint, we are essentially computing the defocusing lengths due to space charge forces here.

      The defocusing constants are computed from a weighted linear regression of the true fields generated by an ellipsoidally symmetric charge distribution. By the equivalent uniform beam principle the space charge effects (to second order) are only loosely coupled to the actual profile of the distribution (assuming that it is ellipsoidally symmetric). The effect from the distribution profile manifests itself as a factor, we take this factor to be that for a uniform ellipsoid for computational purposes.

      Parameters:
      dblGamma - the relativistic factor γ ≡ 1/(1 - v2/c2)1/2 for ellipsoid
      arrMoments - three-array (<x2>, <y2> <z2>) of ellipsoid spatial second moments
      Returns:
      three-array (knx2, kny2, knz2) of defocusing constants
      See Also:
    • compDefocusConstantsAlaTrace3D

      public static double[] compDefocusConstantsAlaTrace3D(double dblGamma, double[] arrMoments)
      This method is provided as a comparison utility for validation against simulation with Trace3D. Trace3D uses an approximation to the elliptic integrals encountered in the space charge field expressions. These approximations break down to first order as the beam ellipsoid becomes increasingly eccentric in the transverse direction. Moreover, there is a preferred direction in these approximations, namely, the axial direction. The approximation is most valid for the situation where the beam is aligned to the axial direction. Thus, for arbitrarily rotated ellipsoids, the approximation becomes less valid.
      Parameters:
      dblGamma - the relativistic factor γ ≡ 1/(1 - v2/c2)1/2 for ellipsoid
      arrMoments - three-array (<xx>,<yy>,<zz>) of ellipsoid spatial second moments
      Returns:
      three-array (knx^2,kny^2,knz^2) of defocusing constants
      See Also:
    • getGamma

      public double getGamma()
      Return the relativistic parameter for the ellipsoidal charge distribution.
      Returns:
      relativistic parameter γ = (1 + v2/c2)1/2
    • getCorrelationLab

      public CovarianceMatrix getCorrelationLab()
      Return the original correlation matrix for the beam in the laboratory coordinates.
      Returns:
      beam correlation matrix in laboratory frame
    • get2ndMomentX

      public double get2ndMomentX()
      Return the value of the first ellipsoid second spatial moment in the stationary beam frame and aligned to the coordinate axes. NOTE The after rotation the (x,y,z) coordinate are somewhat arbitrary and no longer coincide with the original laboratory coordinates.
      Returns:
      second moment <x*x>
    • get2ndMomentY

      public double get2ndMomentY()
      Return the value of the second ellipsoid second spatial moment in the stationary beam frame and aligned to the coordinate axes. NOTE The after rotation the (x,y,z) coordinate are somewhat arbitrary and no longer coincide with the original laboratory coordinates.
      Returns:
      second moment <y*y>
    • get2ndMomentZ

      public double get2ndMomentZ()
      Return the value of the third ellipsoid second spatial moment in the stationary beam frame and aligned to the coordinate axes. NOTE The after rotation the (x,y,z) coordinate are somewhat arbitrary and no longer coincide with the original laboratory coordinates.
      Returns:
      second moment <z*z>
    • get2ndMoments

      public double[] get2ndMoments()
      Return all the ellipsoid second spatial moments in the stationary beam frame and aligned to the coordinate axes. NOTE The after rotation the (x,y,z) coordinate are somewhat arbitrary and no longer coincide with the original laboratory coordinates.
      Returns:
      three-array (<x*x>,<y*y>,<z*z>) second moments
    • getSemiAxisX

      public double getSemiAxisX()
      Return the value of the first ellipsoid semi-axis in the stationary beam frame. The first value of getSemiAxes(). NOTE The after rotation the (x,y,z) coordinate are somewhat arbitrary and no longer coincide with the original laboratory coordinates.
      Returns:
      1st semi-axis value
      See Also:
    • getSemiAxisY

      public double getSemiAxisY()
      Return the value of the second ellipsoid semi-axis in the stationary beam frame. The second value of getSemiAxes(). NOTE The after rotation the (x,y,z) coordinate are somewhat arbitrary and no longer coincide with the original laboratory coordinates.
      Returns:
      2nd semi-axis value
      See Also:
    • getSemiAxisZ

      public double getSemiAxisZ()
      Return the value of the third ellipsoid semi-axis in the stationary beam frame. The third value of getSemiAxes(). NOTE The after rotation the (x,y,z) coordinate are somewhat arbitrary and no longer coincide with the original laboratory coordinates.
      Returns:
      3rd semi-axis value
      See Also:
    • getSemiAxes

      public double[] getSemiAxes()
      Return all the ellipsoid semi-axes lengths as an array. Note that these semi-axes are the values in the stationary beam frame. NOTE The after rotation the (x,y,z) coordinate are somewhat arbitrary and no longer coincide with the original laboratory coordinates.
      Returns:
      array (a,b,c) of ellipsoid semi-axes
    • getDefocusingConstantX

      public double getDefocusingConstantX()
      Return the 1st normalized space charge defocusing constant. NOTE The after rotation the (x,y,z) coordinate are somewhat arbitrary and no longer coincide with the original laboratory coordinates.
      Returns:
      1st normalized defocusing constant knx
      See Also:
    • getDefocusingConstantY

      public double getDefocusingConstantY()
      Return the 2nd normalized space charge defocusing constant. NOTE The after rotation the (x,y,z) coordinate are somewhat arbitrary and no longer coincide with the original laboratory coordinates.
      Returns:
      2nd normalized defocusing constant knx
      See Also:
    • getDefocusingConstantZ

      public double getDefocusingConstantZ()
      Return the 3rd normalized space charge defocusing constant. NOTE The after rotation the (x,y,z) coordinate are somewhat arbitrary and no longer coincide with the original laboratory coordinates.
      Returns:
      3rd normalized defocusing constant knz
      See Also:
    • getDefocusingConstants

      public double[] getDefocusingConstants()

      Return all the normalized space-charge defocusing constants for the beam ellipsoid. These are the inverses of the normalized defocal lengths (fnx,fny,fnz) and are used to construct the space charge generator matrix and space charge transfer matrix. The unnormalized focal lengths f are given by

      f = ds*K*fn

      where ds is the increment path length over which the space charge kick is being applied, K is the generalized (3D) beam perveance, and fn is the normalized defocal length. The normalized defocal length is given by

      fn = 1/kn^2

      where kn^2 is the normalized (squared) focusing constant, i.e., the value returned by this method.

      Returns:
      three-array (knx^2, kny^2, knz^2) of normalized defocusing constants
    • getCorrelationBeam

      public CovarianceMatrix getCorrelationBeam()
      Return the correlation matrix for the beam in the stationary beam coordinates.
      Returns:
      beam correlation matrix in beam frame
    • getLorentzTransform

      public PhaseMatrix getLorentzTransform()
      Get the Lorentz transform matrix which takes the laboratory coordinates to the beam coordinates. Note that the laboratory coordinates actually move with the beam centroid so they are note truly stationary w.r.t. to the machine. Thus, the Lorentz transform provided here does not include any translation effects, only the effects of length contraction and time dilation.
      Returns:
      relativistic Lorentz transform matrix from lab to beam coordinates
    • getTranslation

      public PhaseMatrix getTranslation()
      Return the translation matrix which transforms coordinates in the beam frame to those with the ellipsoid centroid as the coordinate origin.

      Note that because we are using homogeneous coordinates, Galilean translations may be performed with matrix multiplications.
      Returns:
      translation matrix moving the coordinate origin to the beam centroid
    • getRotation

      public PhaseMatrix getRotation()

      Get orthogonal rotation matrix R in SO(7) that rotates the ellipsoid spatial semi-axes onto the spatial coordinate axes.

      The rotation R is actually the Cartesian product of a single rotation r from SO(3) that rotates the ellipsoid's spatial coordinates onto the coordinate axes. That is,

      R = r × r contained in SO(7) contained in R7×7

      In this manner the momentum coordinates at each point (x,y,z) are rotated by an equal amount and so velocities are mapped accordingly.

      Returns:
      rotation matrix in SO(7) moving the ellipsoid into standard position
    • getLabToBeamTransform

      public PhaseMatrix getLabToBeamTransform()
      Return the complete transformation from the laboratory inertial coordinates to the ellipsoid inertial coordinates. The tranform takes coordinates in the laboratory frame to the natural coordinates of the ellipsoid - this coordinate system has the centroid as the origin and the ellipsoid semi-axes are aligned to the coordinate axes. Denoting the returned transformation as M, then it is composed of the following factors

          M = R0*T0*L0

      where L0 is the Lorentz transform into the beam frame, T0 is the Galilean transform to the ellipsoid centroid coordinates, and R0 is the rotation that aligns the ellipsoid semi-axes to the coorinates axes putting it into standard position.
      Returns:
      tranformation taking lab inertial coordinates to beam inertial coordinates
      See Also:
    • getBeamToEllipseTransform

      public PhaseMatrix getBeamToEllipseTransform()
      Return the complete transformation from the beam inertial coordinates to the ellipsoid inertial coordinates. The transform takes coordinates in the beam frame to the natural coordinates of the ellipsoid - this coordinate system has the centroid as the origin and the ellipsoid semi-axes are aligned to the coordinate axes. Denoting the returned transformation as M, then it is composed of the following factors

          M = R0*T0

      where T0 is the Galilean transform to the ellipsoid centroid coordinates, and R0 is the rotation that aligns the ellipsoid semi-axes to the coordinates axes putting it into standard position.
      Returns:
      transformation taking beam inertial coordinates to ellipse inertial coordinates
      See Also:
    • computeScheffGenerator

      public PhaseMatrix computeScheffGenerator(double dblPerveance)

      Calculates the transfer matrix generator for space charge effects from this BeamEllipsoid object given the generalized beam three-dimensional perveance dblPerveance.

      Denoting the returned generator matrix as G then the actual transfer matrix M(s) for the space charge effect is given as

          M(s) = esG

      where s is the incremental path length for the dynamics.

      Note that to obtain this matrix a linear fit to the true fields was performed where the regression is weighted by the distribution itself. According the "Equivalent Beam" principle by Sacherar this regression then only loosely couples to the actual distribution profile of the the ellipsoidal charge. For computational purposes we have assumed a uniform density ellipsoid, but that is of no real consequence practically.

      The resulting generator matrix is the product of the generator matrix G0 in the ellipsoid coordinate (which has a very simple form), and a series of matrix transforms. We have

          G = (R0T0L0)-1 G0 (R0T0L0)

      where L0 is the Lorentz transform into the beam frame, T0 is the Galilean transform to the ellipsoid centroid coordinates, and R0 is the rotation that aligns the ellipsoid semi-axes to the coorinates axes putting it into standard position.

      NOTES:

      One should provide the three-dimensional value for the generalized beam perveance. This value K is defined as

          K = (Q/(2πε0)) *(1/(γ3β2)) *(q/(mc2))

      where Q is the total beam charge, ε0 is the permittivity of free space, γ is the relativistic factor, β is the normalized velocity (to the speed of light), q is the unit charge, m is the beam particle mass, and c is the speed of light.

      Parameters:
      dblPerveance - K, the generalized three-dimensional beam perveance
      Returns:
      G, the transfer matrix generator representing linear space charge effects
      See Also:
    • computeScheffMatrix

      public PhaseMatrix computeScheffMatrix(double dblLen, double dblPerveance)

      Compute and return the transfer matrix for space charge effects due to this beam ellipsoid for the given incremental path length dblLen and generalized beam dblPerveance.

      Note that the returned matrix M has the form

          M(s) = exp(sG)

      where s is the incremental path length for the dynamics, and G is the generator matrix generator for the transform. This generator matrix is returned by the method computeScheffGenerator(double), thus, we could compute the transfer matrix M simply by exponentiating the returned value. However, because the generator matrix in the ellipsoid coordinates, G0, is idempotent (i.e., G0G0 = 0) it is computationally faster to assemble the transfer matrix as the product

          M(s) = esG = (R0T0L0)-1 (I + dsG0) (R0T0L0)

      where L0 is the Lorentz transform into the beam frame, T0 is the Galilean transform to the ellipsoid centroid coordinates, and R0 is the rotation that aligns the ellipsoid semi-axes to the coordinates axes putting it into standard position.

      NOTES:

      One should provide the three-dimensional value for the generalized beam perveance. This value K is defined as

          K = (Q/(2πε0)) *(1/(γ3β2)) *(q/(mc2))

      where Q is the total beam charge, ε0 is the permittivity of free space, γ is the relativistic factor, β is the normalized velocity (to the speed of light), q is the unit charge, m is the beam particle mass, and c is the speed of light.

      Parameters:
      dblLen - ds, the incremental path length for the space charge effects
      dblPerveance - K the generalized three-dimensional beam perveance
      Returns:
      transfer matrix representing linear space charge effects
      See Also:
    • computeDCScheffMatrix

      public PhaseMatrix computeDCScheffMatrix(double dblLen, double dblPerveance)