Skip to content

Particles

beamphysics.ParticleGroup

ParticleGroup(h5=None, data=None)

Particle Group class

Initialized on on openPMD beamphysics particle group:

  • h5: open h5 handle, or str that is a file
  • data: raw data

The required bunch data is stored in .data with keys

  • np.array: x, px, y, py, z, pz, t, status, weight
  • str: species

where:

  • x, y, z are positions in units of [m]
  • px, py, pz are momenta in units of [eV/c]
  • t is time in [s]
  • weight is the macro-charge weight in [C], used for all statistical calulations.
  • species is a proper species name: 'electron', etc.

Optional data:

  • np.array: id

where id is a list of unique integers that identify the particles.

Derived data can be computed as attributes:

  • .gamma, .beta, .beta_x, .beta_y, .beta_z: relativistic factors [1].
  • .r, .theta: cylidrical coordinates [m], [1]
  • .pr, .ptheta: momenta in the radial and angular coordinate directions [eV/c]
  • .Lz: angular momentum about the z axis [m*eV/c]
  • .energy : total energy [eV]
  • .kinetic_energy: total energy - mc^2 in [eV].
  • .higher_order_energy: total energy with quadratic fit in z or t subtracted [eV]
  • .p: total momentum in [eV/c]
  • .mass: rest mass in [eV]
  • .xp, .yp: Slopes \(x' = dx/dz = p_x/p_z\) and \(y' = dy/dz = p_y/p_z\) [1].

Normalized transvere coordinates can also be calculated as attributes:

  • .x_bar, .px_bar, .y_bar, .py_bar in [sqrt(m)]

The normalization is automatically calculated from the covariance matrix. See functions in .statistics for more advanced usage.

Their cooresponding amplitudes are:

.Jx, .Jy [m]

where Jx = (x_bar^2 + px_bar^2 )/2.

The momenta are normalized by the mass, so that <Jx> = norm_emit_x and similar for y.

Statistics of any of these are calculated with:

  • .min(X)
  • .max(X)
  • .ptp(X)
  • .avg(X)
  • .std(X)
  • .cov(X, Y, ...)
  • .histogramdd(X, Y, ..., bins=10, range=None)

with a string X as the name any of the properties above.

Useful beam physics quantities are given as attributes:

  • .norm_emit_x
  • .norm_emit_y
  • .norm_emit_4d
  • .higher_order_energy_spread
  • .average_current

Twiss parameters, including dispersion, for the \(x\) or \(y\) plane:

  • .twiss(plane='x', fraction=0.95, p0C=None)

For convenience, plane='xy' will calculate twiss for both planes.

Twiss matched particles, using a simple linear transformation:

  • .twiss_match(self, beta=None, alpha=None, plane='x', p0c=None, inplace=False)

The weight is required and must sum to > 0. The sum of the weights is in .charge. This can also be set: .charge = 1.234 # C, will rescale the .weight array

All attributes can be accessed with brackets: [key]

Additional keys are allowed for convenience: ['min_prop'] will return .min('prop') ['max_prop'] will return .max('prop') ['ptp_prop'] will return .ptp('prop') ['mean_prop'] will return .avg('prop') ['sigma_prop'] will return .std('prop') ['cov_prop1__prop2'] will return .cov('prop1', 'prop2')[0,1]

Units for all attributes can be accessed by:

  • .units(key)

Particles are often stored at the same time (i.e. from a t-based code), or with the same z position (i.e. from an s-based code.) Routines:

  • drift_to_z(z0)
  • drift_to_t(t0)

help to convert these. If no argument is given, particles will be drifted to the mean. Related properties are:

  • .in_t_coordinates returns True if all particles have the same \(t\) corrdinate
  • .in_z_coordinates returns True if all particles have the same \(z\) corrdinate

Convenient plotting is provided with:

  • .plot(...)
  • .slice_plot(...)

Use help(ParticleGroup.plot), etc. for usage.

Parameters:

Name Type Description Default
h5 str, pathlib.Path, h5py.File, or dict

Source of particle data. Can be a filename (str or Path) to an HDF5 file, an open h5py.File handle, or a dict-like object with openPMD particle data. Mutually exclusive with data.

None
data dict

Raw particle data as a dict with keys 'x', 'px', 'y', 'py', 'z', 'pz', 't', 'status', 'weight', and 'species'. Mutually exclusive with h5.

None

Raises:

Type Description
NotImplementedError

If both h5 and data are provided.

Source code in beamphysics/particles.py
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
def __init__(self, h5=None, data=None):
    if h5 and data:
        # TODO:
        # Allow merging or changing some array with extra data
        raise NotImplementedError("Cannot init on both h5 and data")

    if h5:
        # Allow filename
        if isinstance(h5, (str, pathlib.Path)):
            fname = os.path.expandvars(h5)
            assert os.path.exists(fname), f"File does not exist: {fname}"

            with File(fname, "r") as hh5:
                pp = particle_paths(hh5)
                assert len(pp) == 1, f"Number of particle paths in {h5}: {len(pp)}"
                data = load_bunch_data(hh5[pp[0]])

        elif isinstance(h5, File):
            pp = particle_paths(h5)
            assert len(pp) == 1, f"Number of particle paths in {h5}: {len(pp)}"
            data = load_bunch_data(h5[pp[0]])
        else:
            # Try dict
            data = load_bunch_data(h5)
    else:
        # Fill out data. Exclude species.
        data = full_data(data)
        species = list(set(data["species"]))

        # Allow for empty data (len=0). Otherwise, check species.
        if len(species) >= 1:
            assert len(species) == 1, f"mixed species are not allowed: {species}"
            data["species"] = species[0]

    self._settable_array_keys = [
        "x",
        "px",
        "y",
        "py",
        "z",
        "pz",
        "t",
        "status",
        "weight",
    ]
    # Optional data
    for k in ["id"]:
        if k in data:
            self._settable_array_keys.append(k)

    self._settable_scalar_keys = ["species"]
    self._settable_keys = self._settable_array_keys + self._settable_scalar_keys
    # Internal data. Only allow settable keys
    self._data = {k: data[k] for k in self._settable_keys}

Attributes

beamphysics.ParticleGroup.Jx property
Jx

Normalized amplitude J in the x-px plane

beamphysics.ParticleGroup.Jy property
Jy

Normalized amplitude J in the y-py plane

beamphysics.ParticleGroup.Lz property
Lz

Angular momentum around the z axis in m*eV/c Lz = x * py - y * px

beamphysics.ParticleGroup.average_current property
average_current

Simple average current = charge / dt in [A], with dt = (max_t - min_t) If particles are in \(t\) coordinates, will trydt = (max_z - min_z)*c_light*beta_z

beamphysics.ParticleGroup.beta property
beta

Relativistic beta

beamphysics.ParticleGroup.beta_x property writable
beta_x

Relativistic beta, x component

beamphysics.ParticleGroup.beta_y property writable
beta_y

Relativistic beta, y component

beamphysics.ParticleGroup.beta_z property writable
beta_z

Relativistic beta, z component

beamphysics.ParticleGroup.charge property writable
charge

Total charge in C

beamphysics.ParticleGroup.data property
data

Internal data dict

beamphysics.ParticleGroup.energy property
energy

Total energy in eV

beamphysics.ParticleGroup.gamma property writable
gamma

Relativistic gamma

beamphysics.ParticleGroup.higher_order_energy property
higher_order_energy

Fits a quadratic (order=2) to the Energy vs. time, and returns the energy with this subtracted.

beamphysics.ParticleGroup.higher_order_energy_spread property
higher_order_energy_spread

Legacy syntax to compute the standard deviation of higher_order_energy.

beamphysics.ParticleGroup.id property writable
id

id integer

beamphysics.ParticleGroup.in_t_coordinates property
in_t_coordinates

Returns True if all particles have the same t coordinate

beamphysics.ParticleGroup.in_z_coordinates property
in_z_coordinates

Returns True if all particles have the same z coordinate

beamphysics.ParticleGroup.kinetic_energy property
kinetic_energy

Kinetic energy in eV

beamphysics.ParticleGroup.mass property
mass

Rest mass in eV

beamphysics.ParticleGroup.n_alive property
n_alive

Number of alive particles, defined by status == 1

beamphysics.ParticleGroup.n_dead property
n_dead

Number of alive particles, defined by status != 1

beamphysics.ParticleGroup.n_particle property
n_particle

Total number of particles. Same as len

beamphysics.ParticleGroup.norm_emit_4d property
norm_emit_4d

Normalized emittance in the xy planes (4D)

beamphysics.ParticleGroup.norm_emit_x property
norm_emit_x

Normalized emittance in the x plane

beamphysics.ParticleGroup.norm_emit_y property
norm_emit_y

Normalized emittance in the x plane

beamphysics.ParticleGroup.p property
p

Total momemtum in eV/c

beamphysics.ParticleGroup.pr property
pr

Momentum in the radial direction in eV/c r_hat = cos(theta) xhat + sin(theta) yhat pr = p dot r_hat

beamphysics.ParticleGroup.ptheta property
ptheta

Momentum in the polar theta direction. theta_hat = -sin(theta) xhat + cos(theta) yhat ptheta = p dot theta_hat Note that Lz = r*ptheta

beamphysics.ParticleGroup.px property writable
px

px coordinate in [eV/c]

beamphysics.ParticleGroup.px_bar property
px_bar

Normalized px in units of sqrt(m)

beamphysics.ParticleGroup.py property writable
py

py coordinate in [eV/c]

beamphysics.ParticleGroup.py_bar property
py_bar

Normalized py in units of sqrt(m)

beamphysics.ParticleGroup.pz property writable
pz

pz coordinate in [eV/c]

beamphysics.ParticleGroup.r property
r

Radius in the xy plane: r = sqrt(x^2 + y^2) in m

beamphysics.ParticleGroup.species property writable
species

species string

beamphysics.ParticleGroup.species_charge property
species_charge

Species charge in C

beamphysics.ParticleGroup.status property writable
status

status coordinate in [1]

beamphysics.ParticleGroup.t property writable
t

t coordinate in [s]

beamphysics.ParticleGroup.theta property
theta

Angle in xy plane: theta = arctan2(y, x) in radians

beamphysics.ParticleGroup.weight property writable
weight

weight coordinate in [C]

beamphysics.ParticleGroup.x property writable
x

x coordinate in [m]

beamphysics.ParticleGroup.x_bar property
x_bar

Normalized x in units of sqrt(m)

beamphysics.ParticleGroup.xp property
xp

x slope px/pz (dimensionless)

beamphysics.ParticleGroup.y property writable
y

y coordinate in [m]

beamphysics.ParticleGroup.y_bar property
y_bar

Normalized y in units of sqrt(m)

beamphysics.ParticleGroup.yp property
yp

y slope py/pz (dimensionless)

beamphysics.ParticleGroup.z property writable
z

z coordinate in [m]

Functions

beamphysics.ParticleGroup.apply_wakefield
apply_wakefield(wakefield, length, inplace=False, include_self_kick=True)

Apply wakefield momentum kicks to this ParticleGroup.

Parameters:

Name Type Description Default
wakefield WakefieldBase

A wakefield object providing the particle_kicks(z, weight) method.

required
length float

Length over which the wakefield acts [m].

required
inplace bool

If True, modifies in place. If False, returns a modified copy. Default is False.

False
include_self_kick bool

Whether to include the self-kick term. Default is True.

True

Returns:

Type Description
ParticleGroup or None

Modified ParticleGroup if inplace=False, otherwise None.

Examples:

::

from beamphysics.wakefields import ResistiveWallWakefield
wake = ResistiveWallWakefield.from_material("copper-slac-pub-10707", radius=2.5e-3)
P_after = P.apply_wakefield(wake, length=10.0)
Source code in beamphysics/particles.py
1545
1546
1547
1548
1549
1550
1551
1552
1553
1554
1555
1556
1557
1558
1559
1560
1561
1562
1563
1564
1565
1566
1567
1568
1569
1570
1571
1572
1573
1574
1575
1576
1577
1578
1579
1580
1581
1582
1583
1584
1585
1586
1587
1588
1589
1590
1591
1592
1593
1594
1595
1596
def apply_wakefield(
    self,
    wakefield: WakefieldBase,
    length: float,
    inplace: bool = False,
    include_self_kick: bool = True,
):
    """
    Apply wakefield momentum kicks to this ParticleGroup.

    Parameters
    ----------
    wakefield : WakefieldBase
        A wakefield object providing the `particle_kicks(z, weight)` method.
    length : float
        Length over which the wakefield acts [m].
    inplace : bool, optional
        If True, modifies in place. If False, returns a modified copy.
        Default is False.
    include_self_kick : bool, optional
        Whether to include the self-kick term. Default is True.

    Returns
    -------
    ParticleGroup or None
        Modified ParticleGroup if inplace=False, otherwise None.

    Examples
    --------
    ::

        from beamphysics.wakefields import ResistiveWallWakefield
        wake = ResistiveWallWakefield.from_material("copper-slac-pub-10707", radius=2.5e-3)
        P_after = P.apply_wakefield(wake, length=10.0)
    """
    if not inplace:
        P = self.copy()
    else:
        P = self

    # Extract z positions
    if P.in_t_coordinates:
        z = np.asarray(P.z)
    else:
        z = -c_light * np.asarray(P.t)

    weight = np.asarray(P.weight)
    kicks = wakefield.particle_kicks(z, weight, include_self_kick=include_self_kick)
    P.pz += kicks * length

    if not inplace:
        return P
beamphysics.ParticleGroup.assign_id
assign_id()

Assign unique integer IDs to all particles.

Sets self.id to integers from 1 to n_particle.

Source code in beamphysics/particles.py
383
384
385
386
387
388
389
390
391
def assign_id(self):
    """
    Assign unique integer IDs to all particles.

    Sets ``self.id`` to integers from 1 to ``n_particle``.
    """
    if "id" not in self._settable_array_keys:
        self._settable_array_keys.append("id")
    self.id = np.arange(1, self["n_particle"] + 1)
beamphysics.ParticleGroup.avg
avg(key)

Return the weighted average of a particle property.

Parameters:

Name Type Description Default
key str

Any valid particle property name.

required

Returns:

Type Description
float

Weighted average using self.weight as weights.

Source code in beamphysics/particles.py
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
def avg(self, key: str):
    """
    Return the weighted average of a particle property.

    Parameters
    ----------
    key : str
        Any valid particle property name.

    Returns
    -------
    float
        Weighted average using ``self.weight`` as weights.
    """
    dat = self[key]  # equivalent to self.key for accessing properties above
    if np.isscalar(dat):
        return dat
    return np.average(dat, weights=self.weight)
beamphysics.ParticleGroup.bunching
bunching(wavelength)

Calculate the normalized bunching parameter, which is the magnitude of the complex sum of weighted exponentials at a given point.

The formula for bunching is given by:

\[ B(z, \lambda) = \frac{\left|\sum w_i e^{i k z_i}\right|}{\sum w_i} \]

where: - \( z \) is the position array, - \( \lambda \) is the wavelength, - \( k = \frac{2\pi}{\lambda} \) is the wave number, - \( w_i \) are the weights.

Parameters:

Name Type Description Default
wavelength float

Wavelength of the wave.

required

Returns:

Type Description
complex

The normalized bunching parameter.

Raises:

Type Description
ValueError

If wavelength is not a positive number.

Source code in beamphysics/particles.py
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
def bunching(self, wavelength):
    r"""
    Calculate the normalized bunching parameter, which is the magnitude of the
    complex sum of weighted exponentials at a given point.

    The formula for bunching is given by:

    $$
    B(z, \lambda) = \frac{\left|\sum w_i e^{i k z_i}\right|}{\sum w_i}
    $$

    where:
    - \( z \) is the position array,
    - \( \lambda \) is the wavelength,
    - \( k = \frac{2\pi}{\lambda} \) is the wave number,
    - \( w_i \) are the weights.

    Parameters
    ----------
    wavelength : float
        Wavelength of the wave.


    Returns
    -------
    complex
        The normalized bunching parameter.

    Raises
    ------
    ValueError
        If `wavelength` is not a positive number.
    """

    if self.in_z_coordinates:
        # Approximate z
        z = self.t * self.avg("beta_z") * c_light
    else:
        z = self.z

    return statistics.bunching(z, wavelength, weight=self.weight)
beamphysics.ParticleGroup.copy
copy()

Return a deep copy of this ParticleGroup.

Returns:

Type Description
ParticleGroup

An independent copy with no shared data.

Source code in beamphysics/particles.py
1747
1748
1749
1750
1751
1752
1753
1754
1755
1756
def copy(self):
    """
    Return a deep copy of this ParticleGroup.

    Returns
    -------
    ParticleGroup
        An independent copy with no shared data.
    """
    return deepcopy(self)
beamphysics.ParticleGroup.cov
cov(*keys)

Covariance matrix from any properties

Example: P = ParticleGroup(h5) P.cov('x', 'px', 'y', 'py')

Source code in beamphysics/particles.py
729
730
731
732
733
734
735
736
737
738
739
def cov(self, *keys):
    """
    Covariance matrix from any properties

    Example:
    P = ParticleGroup(h5)
    P.cov('x', 'px', 'y', 'py')

    """
    dats = np.array([self[key] for key in keys])
    return np.cov(dats, aweights=self.weight)
beamphysics.ParticleGroup.delta
delta(key)

Return an attribute array relative to its weighted mean.

Parameters:

Name Type Description Default
key str

Any valid particle property name.

required

Returns:

Type Description
ndarray

Array of self[key] - self.avg(key).

Source code in beamphysics/particles.py
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
def delta(self, key: str):
    """
    Return an attribute array relative to its weighted mean.

    Parameters
    ----------
    key : str
        Any valid particle property name.

    Returns
    -------
    numpy.ndarray
        Array of ``self[key] - self.avg(key)``.
    """
    return self[key] - self.avg(key)
beamphysics.ParticleGroup.drift
drift(delta_t)

Drift all particles for a given time interval.

Updates x, y, z, and t in-place using each particle's velocity.

Parameters:

Name Type Description Default
delta_t float or ndarray

Time interval in seconds. Can be a scalar (same for all particles) or an array with one value per particle.

required
Source code in beamphysics/particles.py
 986
 987
 988
 989
 990
 991
 992
 993
 994
 995
 996
 997
 998
 999
1000
1001
def drift(self, delta_t):
    """
    Drift all particles for a given time interval.

    Updates x, y, z, and t in-place using each particle's velocity.

    Parameters
    ----------
    delta_t : float or numpy.ndarray
        Time interval in seconds. Can be a scalar (same for all particles)
        or an array with one value per particle.
    """
    self.x = self.x + self.beta_x * c_light * delta_t
    self.y = self.y + self.beta_y * c_light * delta_t
    self.z = self.z + self.beta_z * c_light * delta_t
    self.t = self.t + delta_t
beamphysics.ParticleGroup.drift_to_t
drift_to_t(t=None)

Drifts all particles to the same t.

If no z is given, particles will be drifted to the average t

Source code in beamphysics/particles.py
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
def drift_to_t(self, t=None):
    """
    Drifts all particles to the same t.

    If no z is given, particles will be drifted to the average t
    """
    if t is None:
        t = self.avg("t")
    dt = t - self.t
    self.drift(dt)
    # Fix t to be exactly this value
    self.t = np.full(self.n_particle, t)
beamphysics.ParticleGroup.drift_to_z
drift_to_z(z=None)

Drift all particles to the same z position.

Each particle is drifted by the time needed to reach z given its longitudinal velocity. After drifting, all z values are set exactly to z.

Parameters:

Name Type Description Default
z float

Target z position in meters. If None, particles are drifted to the mean z.

None
Source code in beamphysics/particles.py
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
def drift_to_z(self, z=None):
    """
    Drift all particles to the same z position.

    Each particle is drifted by the time needed to reach `z` given
    its longitudinal velocity. After drifting, all z values are set
    exactly to `z`.

    Parameters
    ----------
    z : float, optional
        Target z position in meters. If None, particles are drifted
        to the mean z.
    """
    if z is None:
        z = self.avg("z")
    dt = (z - self.z) / (self.beta_z * c_light)
    self.drift(dt)
    # Fix z to be exactly this value
    self.z = np.full(self.n_particle, z)
beamphysics.ParticleGroup.fractional_split
fractional_split(fractions, key)

Split particles based on a given array key and a list of specified fractions or a single fraction.

Parameters:

Name Type Description Default
fractions float or list of float

A fraction or a list of fractions used for splitting the particles. All values must be between 0 and 1 (exclusive).

required
key str

The attribute of particles to be used for sorting and splitting (e.g., 'z' for longitudinal position).

required

Returns:

Name Type Description
particle_groups list of ParticleGroup

A list of ParticleGroup objects, each representing a subset of particles based on the specified fractions.

Description

This function splits the given group of particles into multiple subsets based on the provided attribute (e.g., position). The splits are determined such that each specified fraction of the total particle weights is separated. The function first sorts the particles by the specified key, computes the cumulative sum of weights, and determines the split values. It then returns a list of ParticleGroup objects representing the split subsets.

Source code in beamphysics/particles.py
1687
1688
1689
1690
1691
1692
1693
1694
1695
1696
1697
1698
1699
1700
1701
1702
1703
1704
1705
1706
1707
1708
1709
1710
1711
1712
1713
1714
1715
1716
1717
1718
1719
1720
1721
1722
1723
1724
1725
1726
1727
1728
1729
1730
1731
1732
1733
1734
1735
1736
1737
1738
1739
1740
1741
1742
1743
1744
1745
def fractional_split(self, fractions: Union[float, int, list], key: str):
    """
    Split particles based on a given array key and a list of specified fractions or a single fraction.

    Parameters
    ----------
    fractions : float or list of float
        A fraction or a list of fractions used for splitting the particles. All values must be between 0 and 1 (exclusive).

    key : str
        The attribute of particles to be used for sorting and splitting (e.g., 'z' for longitudinal position).

    Returns
    -------
    particle_groups : list of ParticleGroup
        A list of ParticleGroup objects, each representing a subset of particles based on the specified fractions.

    Description
    -----------
    This function splits the given group of particles into multiple subsets based on the provided attribute (e.g., position).
    The splits are determined such that each specified fraction of the total particle weights is separated.
    The function first sorts the particles by the specified key, computes the cumulative sum of weights,
    and determines the split values. It then returns a list of ParticleGroup objects representing the split subsets.

    """

    # Ensure fractions is a list
    if isinstance(fractions, (float, int)):
        fractions = [fractions]

    # Validate fraction values
    if any(f <= 0 or f >= 1 for f in fractions):
        raise ValueError("All fraction values must be between 0 and 1 (exclusive)")

    # Sort particles by the specified key
    ixs = np.argsort(self[key])
    sorted_particles = self[ixs]

    # Sorted weights
    ws = sorted_particles.weight
    total_weight = np.sum(ws)
    cw = np.cumsum(ws) / total_weight

    # Use vectorized searchsorted to determine split indices
    fractions = np.array(fractions)
    split_indices = np.searchsorted(cw, fractions, side="right")

    # Create ParticleGroup subsets for each split
    particle_groups = []
    previous_index = 0
    for isplit in split_indices:
        particle_groups.append(sorted_particles[previous_index:isplit])
        previous_index = isplit

    # Add the remaining particles to the last group
    if previous_index < len(sorted_particles):
        particle_groups.append(sorted_particles[previous_index:])

    return particle_groups
beamphysics.ParticleGroup.from_bmad classmethod
from_bmad(bmad_dict)

Convert Bmad particle data as a dict to ParticleGroup data.

See: ParticleGroup.to_bmad or particlegroup_to_bmad

Parameters:

Name Type Description Default
bmad_data

Dict with keys: 'x' 'px' 'y' 'py' 'z' 'pz', 'charge' 'species', 'tref' 'state'

required

Returns:

Type Description
ParticleGroup
Source code in beamphysics/particles.py
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
@classmethod
@functools.wraps(bmad.bmad_to_particlegroup_data)
def from_bmad(cls, bmad_dict):
    """
    Convert Bmad particle data as a dict
    to ParticleGroup data.

    See: ParticleGroup.to_bmad or particlegroup_to_bmad

    Parameters
    ----------
    bmad_data: dict
        Dict with keys:
        'x'
        'px'
        'y'
        'py'
        'z'
        'pz',
        'charge'
        'species',
        'tref'
        'state'

    Returns
    -------
    ParticleGroup
    """
    data = bmad.bmad_to_particlegroup_data(bmad_dict)
    return cls(data=data)
beamphysics.ParticleGroup.higher_order_energy_calc
higher_order_energy_calc(order=2)

Fits a polynmial with order order to the Energy vs. time, , and returns the energy with this subtracted.

Source code in beamphysics/particles.py
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
def higher_order_energy_calc(self, order=2):
    """
    Fits a polynmial with order `order` to the Energy vs. time, , and returns the energy with this subtracted.
    """
    # order=2
    if self.std("z") < 1e-12:
        # must be at a screen. Use t
        t = self.t
    else:
        # All particles at the same time. Use z to calc t
        t = self.z / c_light
    energy = self.energy

    best_fit_coeffs = np.polynomial.polynomial.polyfit(t, energy, order)
    best_fit = np.polynomial.polynomial.polyval(t, best_fit_coeffs)
    return energy - best_fit
beamphysics.ParticleGroup.histogramdd
histogramdd(*keys, bins=10, range=None)

Wrapper for numpy.histogramdd, but accepts property names as keys.

Computes the multidimensional histogram of keys. Internally uses weights.

Example: P.histogramdd('x', 'y', bins=50) Returns: np.array with shape 50x50, edge list

Source code in beamphysics/particles.py
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
def histogramdd(self, *keys, bins=10, range=None):
    """
    Wrapper for numpy.histogramdd, but accepts property names as keys.

    Computes the multidimensional histogram of keys. Internally uses weights.

    Example:
        P.histogramdd('x', 'y', bins=50)
    Returns:
        np.array with shape 50x50, edge list

    """
    H, edges = np.histogramdd(
        np.array([self[k] for k in list(keys)]).T,
        weights=self.weight,
        bins=bins,
        range=range,
    )

    return H, edges
beamphysics.ParticleGroup.info
info(key)

Return information about a specific statistics key.

Parameters:

Name Type Description Default
key str
required

Returns:

Type Description
dict

Information on the specific statistics item, including the following keys: label, mathlabel, units, description, formula, reference, reference_url, category, dtype, and shape

Source code in beamphysics/particles.py
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
def info(self, key: str):
    """
    Return information about a specific statistics key.

    Parameters
    ----------
    key : str

    Returns
    -------
    dict
        Information on the specific statistics item, including the following keys:
        label, mathlabel, units, description, formula, reference,
        reference_url, category, dtype, and shape
    """
    from .standards.statistics import get_all_statistics_by_key

    return get_all_statistics_by_key()[key]
beamphysics.ParticleGroup.linear_point_transform
linear_point_transform(mat3)

Applies a linear transformation to the particle's spatial coordinates and their conjugate momenta.

The spatial coordinates (x, y, z) are transformed using: [x', y', z'] = mat3 @ [x, y, z]

The conjugate momenta (px, py, pz) are transformed using: [px', py', pz'] = inv(mat3.T) @ [px, py, pz]

Parameters:

Name Type Description Default
mat3 array-like of shape (3, 3)

A 3x3 transformation matrix. Must be convertible to a NumPy array. Coordinates must be ordered as (x, y, z).

required

Raises:

Type Description
ValueError

If mat3 is not a 3x3 matrix.

LinAlgError

If the transformation matrix is singular (i.e., not invertible).

Notes

This method modifies the object in-place.

Source code in beamphysics/particles.py
1815
1816
1817
1818
1819
1820
1821
1822
1823
1824
1825
1826
1827
1828
1829
1830
1831
1832
1833
1834
1835
1836
1837
1838
1839
1840
1841
1842
1843
1844
1845
1846
1847
1848
1849
1850
1851
1852
1853
1854
1855
1856
1857
1858
1859
1860
1861
1862
1863
def linear_point_transform(
    self, mat3: Union[np.ndarray, Sequence[Sequence[float]]]
) -> None:
    """
    Applies a linear transformation to the particle's spatial coordinates and
    their conjugate momenta.

    The spatial coordinates (x, y, z) are transformed using:
        [x', y', z'] = mat3 @ [x, y, z]

    The conjugate momenta (px, py, pz) are transformed using:
        [px', py', pz'] = inv(mat3.T) @ [px, py, pz]

    Parameters
    ----------
    mat3 : array-like of shape (3, 3)
        A 3x3 transformation matrix. Must be convertible to a NumPy array.
        Coordinates must be ordered as (x, y, z).

    Raises
    ------
    ValueError
        If `mat3` is not a 3x3 matrix.
    LinAlgError
        If the transformation matrix is singular (i.e., not invertible).

    Notes
    -----
    This method modifies the object in-place.
    """
    mat3 = np.asarray(mat3, dtype=float)
    if mat3.shape != (3, 3):
        raise ValueError("mat3 must be a 3x3 matrix.")

    # Grab original coordinates
    x, y, z = self.x, self.y, self.z
    px, py, pz = self.px, self.py, self.pz

    # Transform spatial coordinates: [x', y', z'] = mat3 @ [x, y, z]
    self.x = mat3[0, 0] * x + mat3[0, 1] * y + mat3[0, 2] * z
    self.y = mat3[1, 0] * x + mat3[1, 1] * y + mat3[1, 2] * z
    self.z = mat3[2, 0] * x + mat3[2, 1] * y + mat3[2, 2] * z

    # Transform momenta using inverse-transpose: [px', py', pz'] = inv(mat3.T) @ [px, py, pz]
    # Calculate inverse-transpose matrix elements
    inv_mat3_T = np.linalg.inv(mat3.T)
    self.px = inv_mat3_T[0, 0] * px + inv_mat3_T[0, 1] * py + inv_mat3_T[0, 2] * pz
    self.py = inv_mat3_T[1, 0] * px + inv_mat3_T[1, 1] * py + inv_mat3_T[1, 2] * pz
    self.pz = inv_mat3_T[2, 0] * px + inv_mat3_T[2, 1] * py + inv_mat3_T[2, 2] * pz
beamphysics.ParticleGroup.max
max(key)

Return the maximum value of a particle property.

Parameters:

Name Type Description Default
key str

Any valid particle property name.

required

Returns:

Type Description
float

Maximum value.

Source code in beamphysics/particles.py
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
def max(self, key: str):
    """
    Return the maximum value of a particle property.

    Parameters
    ----------
    key : str
        Any valid particle property name.

    Returns
    -------
    float
        Maximum value.
    """
    return np.max(self[key])
beamphysics.ParticleGroup.min
min(key)

Return the minimum value of a particle property.

Parameters:

Name Type Description Default
key str

Any valid particle property name.

required

Returns:

Type Description
float

Minimum value.

Source code in beamphysics/particles.py
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
def min(self, key: str):
    """
    Return the minimum value of a particle property.

    Parameters
    ----------
    key : str
        Any valid particle property name.

    Returns
    -------
    float
        Minimum value.
    """
    return np.min(self[key])
beamphysics.ParticleGroup.plot
plot(key1='x', key2=None, bins=None, *, xlim=None, ylim=None, return_figure=False, tex=True, nice=True, ellipse=False, **kwargs)

1d or 2d density plot.

If one key is given, this will plot the density of that key. Example: .plot('x')

If two keys arg given, this will plot a 2d marginal plot. Example: .plot('x', 'px')

Parameters:

Name Type Description Default
particle_group

The object to plot

required
key1

Key to bin on the x-axis

'x'
key2

Key to bin on the y-axis.

None
bins

Number of bins. If None, this will use a heuristic: bins = sqrt(n_particle/4)

None
xlim

Manual setting of the x-axis limits. Note that these are in raw, unscaled units.

None
ylim

Manual setting of the y-axis limits. Note that these are in raw, unscaled units.

None
tex

Use TEX for labels

True
nice

Scale to nice units

True
ellipse

If True, plot an ellipse representing the 2x2 sigma matrix

False
return_figure

If true, return a matplotlib.figure.Figure object

False
**kwargs

Any additional kwargs to send to the the plot in: plt.subplots(**kwargs)

{}

Returns:

Type Description
None or fig: matplotlib.figure.Figure

This only returns a figure object if return_figure=T, otherwise returns None

Source code in beamphysics/particles.py
1350
1351
1352
1353
1354
1355
1356
1357
1358
1359
1360
1361
1362
1363
1364
1365
1366
1367
1368
1369
1370
1371
1372
1373
1374
1375
1376
1377
1378
1379
1380
1381
1382
1383
1384
1385
1386
1387
1388
1389
1390
1391
1392
1393
1394
1395
1396
1397
1398
1399
1400
1401
1402
1403
1404
1405
1406
1407
1408
1409
1410
1411
1412
1413
1414
1415
1416
1417
1418
1419
1420
1421
1422
1423
1424
1425
1426
1427
1428
1429
1430
1431
1432
1433
1434
1435
1436
1437
1438
1439
def plot(
    self,
    key1="x",
    key2=None,
    bins=None,
    *,
    xlim=None,
    ylim=None,
    return_figure=False,
    tex=True,
    nice=True,
    ellipse=False,
    **kwargs,
):
    """
    1d or 2d density plot.

    If one key is given, this will plot the density of that key.
    Example:
        .plot('x')

    If two keys arg given, this will plot a 2d marginal plot.
    Example:
        .plot('x', 'px')


    Parameters
    ----------
    particle_group: ParticleGroup
        The object to plot

    key1: str, default = 't'
        Key to bin on the x-axis

    key2: str, default = None
        Key to bin on the y-axis.

    bins: int, default = None
       Number of bins. If None, this will use a heuristic: bins = sqrt(n_particle/4)

    xlim: tuple, default = None
        Manual setting of the x-axis limits. Note that these are in raw, unscaled units.

    ylim: tuple, default = None
        Manual setting of the y-axis limits. Note that these are in raw, unscaled units.

    tex: bool, default = True
        Use TEX for labels

    nice: bool, default = True
        Scale to nice units

    ellipse: bool, default = True
        If True, plot an ellipse representing the
        2x2 sigma matrix

    return_figure: bool, default = False
        If true, return a matplotlib.figure.Figure object

    **kwargs
        Any additional kwargs to send to the the plot in: plt.subplots(**kwargs)


    Returns
    -------
    None or fig: matplotlib.figure.Figure
        This only returns a figure object if return_figure=T, otherwise returns None

    """

    if not key2:
        fig = density_plot(
            self, key=key1, bins=bins, xlim=xlim, tex=tex, nice=nice, **kwargs
        )
    else:
        fig = marginal_plot(
            self,
            key1=key1,
            key2=key2,
            bins=bins,
            xlim=xlim,
            ylim=ylim,
            tex=tex,
            nice=nice,
            ellipse=ellipse,
            **kwargs,
        )

    if return_figure:
        return fig
beamphysics.ParticleGroup.ptp
ptp(key)

Return the peak-to-peak (max - min) range of a particle property.

Parameters:

Name Type Description Default
key str

Any valid particle property name.

required

Returns:

Type Description
float

Peak-to-peak value.

Source code in beamphysics/particles.py
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
def ptp(self, key: str):
    """
    Return the peak-to-peak (max - min) range of a particle property.

    Parameters
    ----------
    key : str
        Any valid particle property name.

    Returns
    -------
    float
        Peak-to-peak value.
    """
    return np.ptp(self[key])
beamphysics.ParticleGroup.rotate
rotate(*, x_rot=0.0, y_rot=0.0, z_rot=0.0, order='zxy', xc=0.0, yc=0.0, zc=0.0)

Rotate the beam inplace by the specified angles first around the z axis, then the x axis, finally around the y axis.

Parameters:

Name Type Description Default
x_rot float

Rotation around the x axis, radians

0.0
y_rot float

Rotation around the y axis, radians

0.0
z_rot float

Rotation around the z axis, radians

0.0
order str

A 3-character string specifying the rotation order (e.g., 'zxy', 'zyx'). Each character must be one of 'x', 'y', or 'z'.

'zxy'
xc float

Center of rotation, x coordinate. Default to zero

0.0
yc float

Center of rotation, x coordinate. Default to zero

0.0
zc float

Center of rotation, x coordinate. Default to zero

0.0
Notes

This method modifies the object in-place.

Source code in beamphysics/particles.py
1865
1866
1867
1868
1869
1870
1871
1872
1873
1874
1875
1876
1877
1878
1879
1880
1881
1882
1883
1884
1885
1886
1887
1888
1889
1890
1891
1892
1893
1894
1895
1896
1897
1898
1899
1900
1901
1902
1903
1904
1905
1906
1907
1908
1909
1910
1911
1912
1913
1914
1915
1916
1917
1918
1919
1920
1921
1922
def rotate(
    self,
    *,
    x_rot: float = 0.0,
    y_rot: float = 0.0,
    z_rot: float = 0.0,
    order: str = "zxy",
    xc: float = 0.0,
    yc: float = 0.0,
    zc: float = 0.0,
) -> None:
    """
    Rotate the beam inplace by the specified angles first around the z axis, then the x axis, finally around the y axis.

    Parameters
    ----------
    x_rot : float, optional
        Rotation around the x axis, radians
    y_rot : float, optional
        Rotation around the y axis, radians
    z_rot : float, optional
        Rotation around the z axis, radians
    order : str
        A 3-character string specifying the rotation order (e.g., 'zxy', 'zyx').
        Each character must be one of 'x', 'y', or 'z'.
    xc : float
        Center of rotation, x coordinate. Default to zero
    yc : float
        Center of rotation, x coordinate. Default to zero
    zc : float
        Center of rotation, x coordinate. Default to zero

    Notes
    -----
    This method modifies the object in-place.
    """
    # Special case transformations
    if (x_rot == 0.0) and (y_rot == 0.0):
        return self.rotate_z(z_rot, xc=xc, yc=yc)
    if (x_rot == 0.0) and (z_rot == 0.0):
        return self.rotate_y(y_rot, xc=xc, zc=zc)
    if (z_rot == 0.0) and (y_rot == 0.0):
        return self.rotate_x(x_rot, yc=yc, zc=zc)

    # Shift beam
    self.x = self.x - xc
    self.y = self.y - yc
    self.z = self.z - zc

    # Perform rotation
    self.linear_point_transform(
        get_rotation_matrix(x_rot=x_rot, y_rot=y_rot, z_rot=z_rot, order=order)
    )

    # Shift back
    self.x = self.x + xc
    self.y = self.y + yc
    self.z = self.z + zc
beamphysics.ParticleGroup.rotate_x
rotate_x(theta, yc=0.0, zc=0.0)

Rotate the object about the x-axis by angle theta (radians). Affects y–z coordinates and py–pz momenta.

Parameters:

Name Type Description Default
theta float

Angle of rotation, radians

required
yc float

Center of rotation, y coordinate. Default to zero

0.0
zc float

Center of rotation, z coordinate. Default to zero

0.0
Source code in beamphysics/particles.py
1924
1925
1926
1927
1928
1929
1930
1931
1932
1933
1934
1935
1936
1937
1938
1939
1940
1941
1942
1943
1944
1945
1946
1947
1948
1949
1950
1951
1952
1953
1954
1955
1956
1957
1958
1959
1960
1961
def rotate_x(
    self,
    theta: float,
    yc: float = 0.0,
    zc: float = 0.0,
) -> None:
    """
    Rotate the object about the x-axis by angle `theta` (radians).
    Affects y–z coordinates and py–pz momenta.

    Parameters
    ----------
    theta : float
        Angle of rotation, radians
    yc : float
        Center of rotation, y coordinate. Default to zero
    zc : float
        Center of rotation, z coordinate. Default to zero
    """
    # Shift beam
    self.y = self.y - yc
    self.z = self.z - zc

    c = np.cos(theta)
    s = np.sin(theta)

    y, z = self.y, self.z
    py, pz = self.py, self.pz

    self.y = c * y - s * z
    self.z = s * y + c * z

    self.py = c * py - s * pz
    self.pz = s * py + c * pz

    # Shift back
    self.y = self.y + yc
    self.z = self.z + zc
beamphysics.ParticleGroup.rotate_y
rotate_y(theta, xc=0.0, zc=0.0)

Rotate the object about the y-axis by angle theta (radians). Affects z–x coordinates and pz–px momenta.

Parameters:

Name Type Description Default
theta float

Angle of rotation, radians

required
xc float

Center of rotation, x coordinate. Default to zero

0.0
zc float

Center of rotation, z coordinate. Default to zero

0.0
Source code in beamphysics/particles.py
1963
1964
1965
1966
1967
1968
1969
1970
1971
1972
1973
1974
1975
1976
1977
1978
1979
1980
1981
1982
1983
1984
1985
1986
1987
1988
1989
1990
1991
1992
1993
1994
1995
1996
1997
1998
1999
2000
def rotate_y(
    self,
    theta: float,
    xc: float = 0.0,
    zc: float = 0.0,
) -> None:
    """
    Rotate the object about the y-axis by angle `theta` (radians).
    Affects z–x coordinates and pz–px momenta.

    Parameters
    ----------
    theta : float
        Angle of rotation, radians
    xc : float
        Center of rotation, x coordinate. Default to zero
    zc : float
        Center of rotation, z coordinate. Default to zero
    """
    # Shift beam
    self.x = self.x - xc
    self.z = self.z - zc

    c = np.cos(theta)
    s = np.sin(theta)

    z, x = self.z, self.x
    pz, px = self.pz, self.px

    self.z = c * z - s * x
    self.x = s * z + c * x

    self.pz = c * pz - s * px
    self.px = s * pz + c * px

    # Shift beam
    self.x = self.x + xc
    self.z = self.z + zc
beamphysics.ParticleGroup.rotate_z
rotate_z(theta, xc=0.0, yc=0.0)

Rotate the object about the z-axis by angle theta (radians). Affects x–y coordinates and px–py momenta.

Parameters:

Name Type Description Default
theta float

Angle of rotation, radians

required
xc float

Center of rotation, x coordinate. Default to zero

0.0
yc float

Center of rotation, y coordinate. Default to zero

0.0
Source code in beamphysics/particles.py
2002
2003
2004
2005
2006
2007
2008
2009
2010
2011
2012
2013
2014
2015
2016
2017
2018
2019
2020
2021
2022
2023
2024
2025
2026
2027
2028
2029
2030
2031
2032
2033
2034
2035
2036
2037
2038
2039
def rotate_z(
    self,
    theta: float,
    xc: float = 0.0,
    yc: float = 0.0,
) -> None:
    """
    Rotate the object about the z-axis by angle `theta` (radians).
    Affects x–y coordinates and px–py momenta.

    Parameters
    ----------
    theta : float
        Angle of rotation, radians
    xc : float
        Center of rotation, x coordinate. Default to zero
    yc : float
        Center of rotation, y coordinate. Default to zero
    """
    # Shift beam
    self.x = self.x - xc
    self.y = self.y - yc

    c = np.cos(theta)
    s = np.sin(theta)

    x, y = self.x, self.y
    px, py = self.px, self.py

    self.x = c * x - s * y
    self.y = s * x + c * y

    self.px = c * px - s * py
    self.py = s * px + c * py

    # Shift beam
    self.x = self.x + xc
    self.y = self.y + yc
beamphysics.ParticleGroup.slice_plot
slice_plot(*keys, n_slice=100, slice_key=None, tex=True, nice=True, return_figure=False, xlim=None, ylim=None, **kwargs)

Plot slice statistics along the bunch.

Parameters:

Name Type Description Default
*keys str

Keys to plot slice statistics for (e.g., 'sigma_x', 'norm_emit_x').

()
n_slice int

Number of slices. Default is 100.

100
slice_key str

Key to slice on. If None, defaults to 'z' for t-coordinates or 't' for z-coordinates.

None
tex bool

Use TeX for axis labels. Default is True.

True
nice bool

Scale to nice units. Default is True.

True
return_figure bool

If True, return the matplotlib Figure. Default is False.

False
xlim tuple of float

Manual x-axis limits in raw units.

None
ylim tuple of float

Manual y-axis limits in raw units.

None
**kwargs

Additional keyword arguments passed to plt.subplots.

{}

Returns:

Type Description
None or Figure

Returns a Figure only if return_figure=True.

Source code in beamphysics/particles.py
1488
1489
1490
1491
1492
1493
1494
1495
1496
1497
1498
1499
1500
1501
1502
1503
1504
1505
1506
1507
1508
1509
1510
1511
1512
1513
1514
1515
1516
1517
1518
1519
1520
1521
1522
1523
1524
1525
1526
1527
1528
1529
1530
1531
1532
1533
1534
1535
1536
1537
1538
1539
1540
1541
1542
1543
def slice_plot(
    self,
    *keys,
    n_slice=100,
    slice_key=None,
    tex=True,
    nice=True,
    return_figure=False,
    xlim=None,
    ylim=None,
    **kwargs,
):
    """
    Plot slice statistics along the bunch.

    Parameters
    ----------
    *keys : str
        Keys to plot slice statistics for (e.g., ``'sigma_x'``, ``'norm_emit_x'``).
    n_slice : int, optional
        Number of slices. Default is 100.
    slice_key : str, optional
        Key to slice on. If None, defaults to ``'z'`` for t-coordinates or
        ``'t'`` for z-coordinates.
    tex : bool, optional
        Use TeX for axis labels. Default is True.
    nice : bool, optional
        Scale to nice units. Default is True.
    return_figure : bool, optional
        If True, return the matplotlib Figure. Default is False.
    xlim : tuple of float, optional
        Manual x-axis limits in raw units.
    ylim : tuple of float, optional
        Manual y-axis limits in raw units.
    **kwargs
        Additional keyword arguments passed to ``plt.subplots``.

    Returns
    -------
    None or matplotlib.figure.Figure
        Returns a Figure only if ``return_figure=True``.
    """
    fig = slice_plot(
        self,
        *keys,
        n_slice=n_slice,
        slice_key=slice_key,
        tex=tex,
        nice=nice,
        xlim=xlim,
        ylim=ylim,
        **kwargs,
    )

    if return_figure:
        return fig
beamphysics.ParticleGroup.slice_statistics
slice_statistics(*keys, n_slice=100, slice_key=None)

Compute slice statistics along the bunch.

Parameters:

Name Type Description Default
*keys str

Statistical quantities to compute per slice (e.g., 'sigma_x', 'norm_emit_x', 'mean_energy').

()
n_slice int

Number of slices. Default is 100.

100
slice_key str

Key to slice on. If None, defaults to 'z' when particles are in t-coordinates, or 't' otherwise.

None

Returns:

Type Description
dict

Dictionary of slice statistics arrays, including a 'current' or 'density' key computed from charge per slice.

Source code in beamphysics/particles.py
1441
1442
1443
1444
1445
1446
1447
1448
1449
1450
1451
1452
1453
1454
1455
1456
1457
1458
1459
1460
1461
1462
1463
1464
1465
1466
1467
1468
1469
1470
1471
1472
1473
1474
1475
1476
1477
1478
1479
1480
1481
1482
1483
1484
1485
1486
def slice_statistics(self, *keys, n_slice=100, slice_key=None):
    """
    Compute slice statistics along the bunch.

    Parameters
    ----------
    *keys : str
        Statistical quantities to compute per slice (e.g., ``'sigma_x'``,
        ``'norm_emit_x'``, ``'mean_energy'``).
    n_slice : int, optional
        Number of slices. Default is 100.
    slice_key : str, optional
        Key to slice on. If None, defaults to ``'z'`` when particles are in
        t-coordinates, or ``'t'`` otherwise.

    Returns
    -------
    dict
        Dictionary of slice statistics arrays, including a
        ``'current'`` or ``'density'`` key computed from charge
        per slice.
    """

    if slice_key is None:
        if self.in_t_coordinates:
            slice_key = "z"

        else:
            slice_key = "t"

    if slice_key in ("t", "delta_t"):
        density_name = "current"
    else:
        density_name = "density"

    keys = set(keys)
    keys.add("mean_" + slice_key)
    keys.add("ptp_" + slice_key)
    keys.add("charge")
    slice_dat = slice_statistics(
        self, n_slice=n_slice, slice_key=slice_key, keys=keys
    )

    slice_dat[density_name] = slice_dat["charge"] / slice_dat["ptp_" + slice_key]

    return slice_dat
beamphysics.ParticleGroup.split
split(n_chunks=100, key='z')

Split the particle group into even chunks sorted by a key.

Parameters:

Name Type Description Default
n_chunks int

Number of chunks to split into. Default is 100.

100
key str

Particle attribute to sort by before splitting. Default is 'z'.

'z'

Returns:

Type Description
list of ParticleGroup

A list of ParticleGroup objects, one per chunk.

Source code in beamphysics/particles.py
1669
1670
1671
1672
1673
1674
1675
1676
1677
1678
1679
1680
1681
1682
1683
1684
1685
def split(self, n_chunks=100, key="z"):
    """
    Split the particle group into even chunks sorted by a key.

    Parameters
    ----------
    n_chunks : int, optional
        Number of chunks to split into. Default is 100.
    key : str, optional
        Particle attribute to sort by before splitting. Default is ``'z'``.

    Returns
    -------
    list of ParticleGroup
        A list of ParticleGroup objects, one per chunk.
    """
    return split_particles(self, n_chunks=n_chunks, key=key)
beamphysics.ParticleGroup.std
std(key)

Return the weighted standard deviation of a particle property.

Parameters:

Name Type Description Default
key str

Any valid particle property name.

required

Returns:

Type Description
float

Weighted standard deviation using self.weight as weights.

Source code in beamphysics/particles.py
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
def std(self, key: str):
    """
    Return the weighted standard deviation of a particle property.

    Parameters
    ----------
    key : str
        Any valid particle property name.

    Returns
    -------
    float
        Weighted standard deviation using ``self.weight`` as weights.
    """
    dat = self[key]
    if np.isscalar(dat):
        return 0
    avg_dat = self.avg(key)
    return np.sqrt(np.average((dat - avg_dat) ** 2, weights=self.weight))
beamphysics.ParticleGroup.twiss
twiss(plane='x', fraction=1, p0c=None)

Returns Twiss and Dispersion dict.

plane can be:

'x', 'y', or 'xy'

Optionally a fraction of the particles, based on amplitiude, can be specified.

Source code in beamphysics/particles.py
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
def twiss(self, plane="x", fraction=1, p0c=None):
    """
    Returns Twiss and Dispersion dict.

    plane can be:

    `'x'`, `'y'`, or `'xy'`

    Optionally a fraction of the particles, based on amplitiude, can be specified.
    """
    d = {}
    for p in plane:
        d.update(
            particle_twiss_dispersion(self, plane=p, fraction=fraction, p0c=p0c)
        )
    return d
beamphysics.ParticleGroup.twiss_match
twiss_match(beta=None, alpha=None, plane='x', p0c=None, inplace=False)

Returns a ParticleGroup with requested Twiss parameters.

See: statistics.matched_particles

Source code in beamphysics/particles.py
795
796
797
798
799
800
801
802
803
804
def twiss_match(self, beta=None, alpha=None, plane="x", p0c=None, inplace=False):
    """
    Returns a ParticleGroup with requested Twiss parameters.

    See: statistics.matched_particles
    """

    return matched_particles(
        self, beta=beta, alpha=alpha, plane=plane, inplace=inplace
    )
beamphysics.ParticleGroup.units
units(key)

Get the unit information (pmd_unit) for a given key.

Parameters:

Name Type Description Default
key str

Any valid particle key name.

required

Returns:

Type Description
pmd_unit
Source code in beamphysics/particles.py
408
409
410
411
412
413
414
415
416
417
418
419
420
421
def units(self, key: str) -> pmd_unit:
    """
    Get the unit information (`pmd_unit`) for a given key.

    Parameters
    ----------
    key : str
        Any valid particle key name.

    Returns
    -------
    pmd_unit
    """
    return pg_units(key)
beamphysics.ParticleGroup.wakefield_plot
wakefield_plot(wake, key=None, nice=True, ax=None, xlim=None, ylim=None, tex=True, bins=None, **kwargs)

Plot per-particle wakefield kicks overlaid with the bunch density.

This function overlays the computed wakefield kicks (in eV/m) as a scatter plot on the primary y-axis, and the corresponding particle density as a histogram on a secondary y-axis. The independent variable is chosen automatically based on the coordinate system or specified explicitly with key.

Parameters:

Name Type Description Default
wake WakefieldBase

A wakefield object providing the particle_kicks(z, weight) method, returning longitudinal wakefield kicks in eV/m.

required
key str

Key to use as the independent variable. If None, defaults to 'delta_z/c' or 'delta_t' depending on particle_group.in_t_coordinates.

None
nice bool

If True, applies unit-aware scaling using SI prefixes (e.g., mm, ns).

True
ax Axes

An existing Axes to plot into. If None, a new figure and axes are created.

None
xlim tuple of float

Limits to apply to the x-axis, in native units.

None
ylim tuple of float

Limits to apply to the y-axis (wakefield kick), in native units.

None
tex bool

Whether to use TeX-style math formatting for labels.

True
bins int or str

Number of bins to use for the density histogram.

None
kwargs dict

Additional keyword arguments passed to plt.subplots() if a new axis is created.

{}

Returns:

Name Type Description
fig Figure

The matplotlib figure containing the plot.

Source code in beamphysics/particles.py
1598
1599
1600
1601
1602
1603
1604
1605
1606
1607
1608
1609
1610
1611
1612
1613
1614
1615
1616
1617
1618
1619
1620
1621
1622
1623
1624
1625
1626
1627
1628
1629
1630
1631
1632
1633
1634
1635
1636
1637
1638
1639
1640
1641
1642
1643
1644
1645
1646
1647
1648
1649
1650
1651
1652
1653
1654
1655
1656
1657
1658
1659
1660
1661
1662
1663
1664
1665
1666
def wakefield_plot(
    self,
    wake: WakefieldBase,
    key=None,
    nice=True,
    ax=None,
    xlim=None,
    ylim=None,
    tex=True,
    bins=None,
    **kwargs,
):
    """
    Plot per-particle wakefield kicks overlaid with the bunch density.

    This function overlays the computed wakefield kicks (in eV/m) as a scatter plot on
    the primary y-axis, and the corresponding particle density as a histogram on a
    secondary y-axis. The independent variable is chosen automatically based on the
    coordinate system or specified explicitly with `key`.

    Parameters
    ----------
    wake : WakefieldBase
        A wakefield object providing the `particle_kicks(z, weight)` method,
        returning longitudinal wakefield kicks in eV/m.

    key : str, optional
        Key to use as the independent variable. If None, defaults to 'delta_z/c' or 'delta_t'
        depending on `particle_group.in_t_coordinates`.

    nice : bool, default=True
        If True, applies unit-aware scaling using SI prefixes (e.g., mm, ns).

    ax : matplotlib.axes.Axes, optional
        An existing Axes to plot into. If None, a new figure and axes are created.

    xlim : tuple of float, optional
        Limits to apply to the x-axis, in native units.

    ylim : tuple of float, optional
        Limits to apply to the y-axis (wakefield kick), in native units.

    tex : bool, default=True
        Whether to use TeX-style math formatting for labels.

    bins : int or str, optional
        Number of bins to use for the density histogram.

    kwargs : dict
        Additional keyword arguments passed to `plt.subplots()` if a new axis is created.

    Returns
    -------
    fig : matplotlib.figure.Figure
        The matplotlib figure containing the plot.
    """

    wakefield_plot(
        self,
        wake,
        key=key,
        nice=nice,
        ax=ax,
        xlim=xlim,
        ylim=ylim,
        tex=tex,
        bins=bins,
        **kwargs,
    )
beamphysics.ParticleGroup.where
where(x)

Filter particles by a boolean condition.

Parameters:

Name Type Description Default
x array_like of bool

Boolean mask or condition array. Particles where x is True are selected.

required

Returns:

Type Description
ParticleGroup

A new ParticleGroup containing only the selected particles.

Source code in beamphysics/particles.py
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
def where(self, x):
    """
    Filter particles by a boolean condition.

    Parameters
    ----------
    x : array_like of bool
        Boolean mask or condition array. Particles where `x` is True
        are selected.

    Returns
    -------
    ParticleGroup
        A new ParticleGroup containing only the selected particles.
    """
    return self[np.where(x)]
beamphysics.ParticleGroup.write
write(h5, name=None)

Write particle data to an HDF5 file or group in openPMD format.

Parameters:

Name Type Description Default
h5 str or Path or File or Group

Target for writing: - str or pathlib.Path: path to a new HDF5 file to be created (opened in mode "w"). Environment variables in path strings are expanded (via os.path.expandvars) before file creation. - h5py.File: an already opened file handle; a new "particles" group is created. - h5py.Group: an existing group; assumed to be an appropriate location to write data. This is for advanced users who know what they are doing. Any other object is treated as a group-like handle compatible with write_pmd_bunch.

required
name str

Name for the subgroup/bunch written inside the "particles" group (or provided group). If None, write_pmd_bunch will write directly to the "particles" group.

None

Returns:

Type Description
None

Examples:

Write to a new file:

>>> P.write("beam.h5")

Write to an already opened file:

>>> import h5py
>>> with h5py.File("beam.h5", "w") as f:
...     P.write(f)

Write into an existing group and name the bunch:

>>> import h5py
>>> with h5py.File("beam.h5", "w") as f:
...     grp = f["particles"]
...     particles.write(grp, name="bunch1")
Source code in beamphysics/particles.py
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
1307
1308
1309
1310
1311
1312
1313
1314
1315
1316
1317
1318
1319
1320
1321
1322
1323
1324
1325
1326
1327
1328
1329
1330
1331
1332
1333
1334
1335
1336
1337
1338
1339
1340
1341
1342
1343
1344
1345
1346
def write(self, h5, name=None) -> None:
    """
    Write particle data to an HDF5 file or group in openPMD format.

    Parameters
    ----------
    h5 : str or pathlib.Path or h5py.File or h5py.Group
        Target for writing:
        - str or pathlib.Path: path to a new HDF5 file to be created (opened in mode "w").
          Environment variables in path strings are expanded (via
          os.path.expandvars) before file creation.
        - h5py.File: an already opened file handle; a new "particles" group is created.
        - h5py.Group: an existing group; assumed to be an appropriate location to write data.
            This is for advanced users who know what they are doing.
            Any other object is treated as a group-like handle compatible with `write_pmd_bunch`.
    name : str, optional
        Name for the subgroup/bunch written inside the "particles" group (or provided group).
        If None, `write_pmd_bunch` will write directly to the "particles" group.

    Returns
    -------
    None

    Examples
    --------
    Write to a new file:
    >>> P.write("beam.h5")

    Write to an already opened file:
    >>> import h5py
    >>> with h5py.File("beam.h5", "w") as f:
    ...     P.write(f)

    Write into an existing group and name the bunch:
    >>> import h5py
    >>> with h5py.File("beam.h5", "w") as f:
    ...     grp = f["particles"]
    ...     particles.write(grp, name="bunch1")
    """
    if isinstance(h5, (str, pathlib.Path)):
        fname = os.path.expandvars(h5)
        g = File(fname, "w")
        pmd_init(g, basePath="/", particlesPath="particles")
        g = g.create_group("particles")
    elif isinstance(h5, File):
        pmd_init(h5, basePath="/", particlesPath="particles")
        g = h5.create_group("particles")
    else:
        g = h5

    write_pmd_bunch(g, self, name=name)
beamphysics.ParticleGroup.write_bmad
write_bmad(filePath, p0c=None, t_ref=0, verbose=False)

Write particle data in Bmad format.

Parameters:

Name Type Description Default
filePath str or Path

Output file path.

required
p0c float

Reference momentum in eV/c. If None, the average pz is used.

None
t_ref float

Reference time in seconds. Default is 0.

0
verbose bool

If True, print information during writing. Default is False.

False
Source code in beamphysics/particles.py
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
def write_bmad(self, filePath, p0c=None, t_ref=0, verbose=False):
    """
    Write particle data in Bmad format.

    Parameters
    ----------
    filePath : str or pathlib.Path
        Output file path.
    p0c : float, optional
        Reference momentum in eV/c. If None, the average pz is used.
    t_ref : float, optional
        Reference time in seconds. Default is 0.
    verbose : bool, optional
        If True, print information during writing. Default is False.
    """
    bmad.write_bmad(self, filePath, p0c=p0c, t_ref=t_ref, verbose=verbose)
beamphysics.ParticleGroup.write_elegant
write_elegant(filePath, verbose=False)

Write particle data in Elegant SDDS format.

Parameters:

Name Type Description Default
filePath str or Path

Output file path.

required
verbose bool

If True, print information during writing. Default is False.

False
Source code in beamphysics/particles.py
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
def write_elegant(self, filePath, verbose=False):
    """
    Write particle data in Elegant SDDS format.

    Parameters
    ----------
    filePath : str or pathlib.Path
        Output file path.
    verbose : bool, optional
        If True, print information during writing. Default is False.
    """
    write_elegant(self, filePath, verbose=verbose)
beamphysics.ParticleGroup.write_genesis2_beam_file
write_genesis2_beam_file(filePath, n_slice=None, verbose=False)

Write particle data as a Genesis 1.3 v2 beam file.

Parameters:

Name Type Description Default
filePath str or Path

Output file path.

required
n_slice int

Number of slices. If None, a default is determined automatically.

None
verbose bool

If True, print information during writing. Default is False.

False
Source code in beamphysics/particles.py
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
def write_genesis2_beam_file(self, filePath, n_slice=None, verbose=False):
    """
    Write particle data as a Genesis 1.3 v2 beam file.

    Parameters
    ----------
    filePath : str or pathlib.Path
        Output file path.
    n_slice : int, optional
        Number of slices. If None, a default is determined automatically.
    verbose : bool, optional
        If True, print information during writing. Default is False.
    """
    # Get beam columns
    beam_columns = genesis2_beam_data(self, n_slice=n_slice)
    # Actually write the file
    write_genesis2_beam_file(filePath, beam_columns, verbose=verbose)
beamphysics.ParticleGroup.write_genesis4_beam
write_genesis4_beam(filePath, n_slice=None, return_input_str=False, verbose=False)

Write particle data as a Genesis 1.3 v4 beam file.

Parameters:

Name Type Description Default
filePath str or Path

Output file path.

required
n_slice int

Number of slices. If None, a default is determined automatically.

None
return_input_str bool

Defaults to False.

False
verbose bool

If True, print information during writing. Default is False.

False
Source code in beamphysics/particles.py
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
@functools.wraps(write_genesis4_beam)
def write_genesis4_beam(
    self, filePath, n_slice=None, return_input_str=False, verbose=False
):
    """
    Write particle data as a Genesis 1.3 v4 beam file.

    Parameters
    ----------
    filePath : str or pathlib.Path
        Output file path.
    n_slice : int, optional
        Number of slices. If None, a default is determined automatically.
    return_input_str : bool, optional
        Defaults to False.
    verbose : bool, optional
        If True, print information during writing. Default is False.
    """
    return write_genesis4_beam(
        self,
        filePath,
        n_slice=n_slice,
        return_input_str=return_input_str,
        verbose=verbose,
    )
beamphysics.ParticleGroup.write_genesis4_distribution
write_genesis4_distribution(filePath, verbose=False)

Write particle data as a Genesis 1.3 v4 distribution HDF5 file.

Parameters:

Name Type Description Default
filePath str or Path

Output file path.

required
verbose bool

If True, print information during writing. Default is False.

False
Source code in beamphysics/particles.py
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
def write_genesis4_distribution(self, filePath, verbose=False):
    """
    Write particle data as a Genesis 1.3 v4 distribution HDF5 file.

    Parameters
    ----------
    filePath : str or pathlib.Path
        Output file path.
    verbose : bool, optional
        If True, print information during writing. Default is False.
    """
    write_genesis4_distribution(self, filePath, verbose=verbose)
beamphysics.ParticleGroup.write_gpt
write_gpt(filePath, asci2gdf_bin=None, verbose=False)

Write particle data in GPT (General Particle Tracer) format.

Parameters:

Name Type Description Default
filePath str or Path

Output file path.

required
asci2gdf_bin str

Path to the asci2gdf binary for GDF conversion. If None, an ASCII file is written.

None
verbose bool

If True, print information during writing. Default is False.

False
Source code in beamphysics/particles.py
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
def write_gpt(self, filePath, asci2gdf_bin=None, verbose=False):
    """
    Write particle data in GPT (General Particle Tracer) format.

    Parameters
    ----------
    filePath : str or pathlib.Path
        Output file path.
    asci2gdf_bin : str, optional
        Path to the ``asci2gdf`` binary for GDF conversion. If None,
        an ASCII file is written.
    verbose : bool, optional
        If True, print information during writing. Default is False.
    """
    write_gpt(self, filePath, asci2gdf_bin=asci2gdf_bin, verbose=verbose)
beamphysics.ParticleGroup.write_impact
write_impact(filePath, cathode_kinetic_energy_ref=None, include_header=True, verbose=False)

Write particle data in Impact-T format.

Parameters:

Name Type Description Default
filePath str or Path

Output file path.

required
cathode_kinetic_energy_ref float

Reference kinetic energy at the cathode in eV. If None, the average kinetic energy is used.

None
include_header bool

If True, include the Impact-T header line. Default is True.

True
verbose bool

If True, print information during writing. Default is False.

False
Source code in beamphysics/particles.py
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
def write_impact(
    self,
    filePath,
    cathode_kinetic_energy_ref=None,
    include_header=True,
    verbose=False,
):
    """
    Write particle data in Impact-T format.

    Parameters
    ----------
    filePath : str or pathlib.Path
        Output file path.
    cathode_kinetic_energy_ref : float, optional
        Reference kinetic energy at the cathode in eV. If None, the
        average kinetic energy is used.
    include_header : bool, optional
        If True, include the Impact-T header line. Default is True.
    verbose : bool, optional
        If True, print information during writing. Default is False.
    """
    from .interfaces.impact import write_impact

    return write_impact(
        self,
        filePath,
        cathode_kinetic_energy_ref=cathode_kinetic_energy_ref,
        include_header=include_header,
        verbose=verbose,
    )
beamphysics.ParticleGroup.write_litrack
write_litrack(filePath, p0c=None, verbose=False)

Write particle data in LiTrack format.

Parameters:

Name Type Description Default
filePath str or Path

Output file path.

required
p0c float

Reference momentum in eV/c. If None, the average pz is used.

None
verbose bool

If True, print information during writing. Default is False.

False
Source code in beamphysics/particles.py
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
def write_litrack(self, filePath, p0c=None, verbose=False):
    """
    Write particle data in LiTrack format.

    Parameters
    ----------
    filePath : str or pathlib.Path
        Output file path.
    p0c : float, optional
        Reference momentum in eV/c. If None, the average pz is used.
    verbose : bool, optional
        If True, print information during writing. Default is False.
    """
    return write_litrack(self, outfile=filePath, p0c=p0c, verbose=verbose)
beamphysics.ParticleGroup.write_lucretia
write_lucretia(filePath, ele_name='BEGINNING', t_ref=0, stop_ix=None, verbose=False)

Write particle data in Lucretia format.

Parameters:

Name Type Description Default
filePath str or Path

Output file path.

required
ele_name str

Element name for the Lucretia beam. Default is 'BEGINNING'.

'BEGINNING'
t_ref float

Reference time in seconds. Default is 0.

0
stop_ix int

Stop index for the Lucretia beam structure.

None
verbose bool

If True, print information during writing. Default is False.

False
Source code in beamphysics/particles.py
1236
1237
1238
1239
1240
1241
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
def write_lucretia(
    self, filePath, ele_name="BEGINNING", t_ref=0, stop_ix=None, verbose=False
):
    """
    Write particle data in Lucretia format.

    Parameters
    ----------
    filePath : str or pathlib.Path
        Output file path.
    ele_name : str, optional
        Element name for the Lucretia beam. Default is ``'BEGINNING'``.
    t_ref : float, optional
        Reference time in seconds. Default is 0.
    stop_ix : int, optional
        Stop index for the Lucretia beam structure.
    verbose : bool, optional
        If True, print information during writing. Default is False.
    """
    return write_lucretia(
        self, filePath, ele_name=ele_name, t_ref=t_ref, stop_ix=stop_ix
    )
beamphysics.ParticleGroup.write_opal
write_opal(filePath, verbose=False, dist_type='emitted')

Write particle data in OPAL (Object Oriented Parallel Accelerator Library) format.

Parameters:

Name Type Description Default
filePath str or Path

Output file path.

required
verbose bool

If True, print information during writing. Default is False.

False
dist_type str

Distribution type, either 'emitted' or 'injected'. Default is 'emitted'.

'emitted'
Source code in beamphysics/particles.py
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
1289
1290
1291
1292
1293
def write_opal(self, filePath, verbose=False, dist_type="emitted"):
    """
    Write particle data in OPAL (Object Oriented Parallel Accelerator Library) format.

    Parameters
    ----------
    filePath : str or pathlib.Path
        Output file path.
    verbose : bool, optional
        If True, print information during writing. Default is False.
    dist_type : str, optional
        Distribution type, either ``'emitted'`` or ``'injected'``.
        Default is ``'emitted'``.
    """
    return write_opal(self, filePath, verbose=verbose, dist_type=dist_type)
beamphysics.ParticleGroup.write_simion
write_simion(filePath, color=0, flip_z_to_x=True, verbose=False)

Write particle data in SIMION format.

Parameters:

Name Type Description Default
filePath str or Path

Output file path.

required
color int

SIMION color index for particles. Default is 0.

0
flip_z_to_x bool

If True, swap z and x axes for SIMION's coordinate convention. Default is True.

True
verbose bool

If True, print information during writing. Default is False.

False
Source code in beamphysics/particles.py
1259
1260
1261
1262
1263
1264
1265
1266
1267
1268
1269
1270
1271
1272
1273
1274
1275
1276
1277
def write_simion(self, filePath, color=0, flip_z_to_x=True, verbose=False):
    """
    Write particle data in SIMION format.

    Parameters
    ----------
    filePath : str or pathlib.Path
        Output file path.
    color : int, optional
        SIMION color index for particles. Default is 0.
    flip_z_to_x : bool, optional
        If True, swap z and x axes for SIMION's coordinate convention.
        Default is True.
    verbose : bool, optional
        If True, print information during writing. Default is False.
    """
    return write_simion(
        self, filePath, verbose=verbose, color=color, flip_z_to_x=flip_z_to_x
    )