Particles
ParticleGroup
ParticleGroup(h5=None, data=None)
Particle Group class
Initialized on on openPMD beamphysics particle group:
- h5: open h5 handle, or
strthat 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,weightstr:species
where:
x,y,zare positions in units of [m]px,py,pzare momenta in units of [eV/c]tis time in [s]weightis the macro-charge weight in [C], used for all statistical calulations.speciesis 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_barin [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_coordinatesreturnsTrueif all particles have the same \(t\) corrdinate.in_z_coordinatesreturnsTrueif all particles have the same \(z\) corrdinate
Convenient plotting is provided with:
.plot(...).slice_plot(...)
Use help(ParticleGroup.plot), etc. for usage.
Parameters:
-
h5(str, pathlib.Path, h5py.File, or dict, default:None) –Source of particle data. Can be a filename (str or Path) to an HDF5 file, an open
h5py.Filehandle, or a dict-like object with openPMD particle data. Mutually exclusive withdata. -
data(dict, default:None) –Raw particle data as a dict with keys
'x','px','y','py','z','pz','t','status','weight', and'species'. Mutually exclusive withh5.
Raises:
-
NotImplementedError–If both
h5anddataare provided.
Methods:
-
__add__–Overloads the + operator to join particle groups.
-
__contains__–Checks internal data
-
__eq__–Check equality of internal data
-
__getitem__–Returns a property or statistical quantity that can be computed:
-
__len__–Return the number of particles.
-
__repr__–Return a detailed string representation with memory location.
-
__str__–Return a human-readable summary string.
-
apply_wakefield–Apply wakefield momentum kicks to this ParticleGroup.
-
assign_id–Assign unique integer IDs to all particles.
-
avg–Return the weighted average of a particle property.
-
bunching–Calculate the normalized bunching parameter, which is the magnitude of the
-
copy–Return a deep copy of this ParticleGroup.
-
cov–Covariance matrix from any properties
-
delta–Return an attribute array relative to its weighted mean.
-
drift–Drift all particles for a given time interval.
-
drift_to_t–Drifts all particles to the same t.
-
drift_to_z–Drift all particles to the same z position.
-
fractional_split–Split particles based on a given array key and a list of specified fractions or a single fraction.
-
from_bmad–Convert Bmad particle data to a ParticleGroup.
-
higher_order_energy_calc–Fits a polynmial with order
orderto the Energy vs. time, , and returns the energy with this subtracted. -
histogramdd–Wrapper for numpy.histogramdd, but accepts property names as keys.
-
info–Return information about a specific statistics key.
-
linear_point_transform–Applies a linear transformation to the particle's spatial coordinates and
-
max–Return the maximum value of a particle property.
-
min–Return the minimum value of a particle property.
-
plot–1d or 2d density plot.
-
ptp–Return the peak-to-peak (max - min) range of a particle property.
-
resample–Resample particles randomly.
-
rotate–Rotate the beam inplace by the specified angles first around the z axis, then the x axis, finally around the y axis.
-
rotate_x–Rotate the object about the x-axis by angle
theta(radians). -
rotate_y–Rotate the object about the y-axis by angle
theta(radians). -
rotate_z–Rotate the object about the z-axis by angle
theta(radians). -
slice_plot–Plot slice statistics along the bunch.
-
slice_statistics–Compute slice statistics along the bunch.
-
split–Split the particle group into even chunks sorted by a key.
-
std–Return the weighted standard deviation of a particle property.
-
to_bmad–Convert to Bmad phase space coordinates.
-
twiss–Returns Twiss and Dispersion dict.
-
twiss_match–Returns a ParticleGroup with requested Twiss parameters.
-
units–Get the unit information (
pmd_unit) for a given key. -
wakefield_plot–Plot per-particle wakefield kicks overlaid with the bunch density.
-
where–Filter particles by a boolean condition.
-
write–Write particle data to an HDF5 file or group in openPMD format.
-
write_astra–Write Astra style particles.
-
write_bmad–Write particle data in Bmad format.
-
write_elegant–Write particle data in Elegant SDDS format.
-
write_genesis2_beam_file–Write particle data as a Genesis 1.3 v2 beam file.
-
write_genesis4_beam–Write particle data as a Genesis 1.3 v4 beam file.
-
write_genesis4_distribution–Write particle data as a Genesis 1.3 v4 distribution HDF5 file.
-
write_gpt–Write particle data in GPT (General Particle Tracer) format.
-
write_impact–Write particle data in Impact-T format.
-
write_litrack–Write particle data in LiTrack format.
-
write_lucretia–Write particle data in Lucretia format.
-
write_opal–Write particle data in OPAL (Object Oriented Parallel Accelerator Library) format.
-
write_simion–Write particle data in SIMION format.
Attributes:
-
Jx–Normalized amplitude J in the x-px plane
-
Jy–Normalized amplitude J in the y-py plane
-
Lz–Angular momentum around the z axis in m*eV/c
-
average_current–Simple average
current = charge / dtin [A], withdt = (max_t - min_t) -
beta–Relativistic beta
-
beta_x–Relativistic beta, x component
-
beta_y–Relativistic beta, y component
-
beta_z–Relativistic beta, z component
-
charge–Total charge in C
-
data–Internal data dict
-
energy–Total energy in eV
-
gamma–Relativistic gamma
-
higher_order_energy–Fits a quadratic (order=2) to the Energy vs. time, and returns the energy with this subtracted.
-
higher_order_energy_spread–Legacy syntax to compute the standard deviation of higher_order_energy.
-
id–id integer
-
in_t_coordinates–Returns True if all particles have the same t coordinate
-
in_z_coordinates–Returns True if all particles have the same z coordinate
-
kinetic_energy–Kinetic energy in eV
-
mass–Rest mass in eV
-
n_alive–Number of alive particles, defined by status == 1
-
n_dead–Number of alive particles, defined by status != 1
-
n_particle–Total number of particles. Same as len
-
norm_emit_4d–Normalized emittance in the xy planes (4D)
-
norm_emit_x–Normalized emittance in the x plane
-
norm_emit_y–Normalized emittance in the x plane
-
p–Total momemtum in eV/c
-
pr–Momentum in the radial direction in eV/c
-
ptheta–Momentum in the polar theta direction.
-
px–px coordinate in [eV/c]
-
px_bar–Normalized px in units of sqrt(m)
-
py–py coordinate in [eV/c]
-
py_bar–Normalized py in units of sqrt(m)
-
pz–pz coordinate in [eV/c]
-
r–Radius in the xy plane: r = sqrt(x^2 + y^2) in m
-
species–species string
-
species_charge–Species charge in C
-
status–status coordinate in [1]
-
t–t coordinate in [s]
-
theta–Angle in xy plane: theta = arctan2(y, x) in radians
-
weight–weight coordinate in [C]
-
x–x coordinate in [m]
-
x_bar–Normalized x in units of sqrt(m)
-
xp–x slope px/pz (dimensionless)
-
y–y coordinate in [m]
-
y_bar–Normalized y in units of sqrt(m)
-
yp–y slope py/pz (dimensionless)
-
z–z coordinate in [m]
Source code in beamphysics/particles.py
189 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 | |
Jx
property
Jx
Normalized amplitude J in the x-px plane
Jy
property
Jy
Normalized amplitude J in the y-py plane
Lz
property
Lz
Angular momentum around the z axis in m*eV/c Lz = x * py - y * px
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
beta
property
beta
Relativistic beta
beta_x
property
writable
beta_x
Relativistic beta, x component
beta_y
property
writable
beta_y
Relativistic beta, y component
beta_z
property
writable
beta_z
Relativistic beta, z component
charge
property
writable
charge
Total charge in C
data
property
data
Internal data dict
energy
property
energy
Total energy in eV
gamma
property
writable
gamma
Relativistic gamma
higher_order_energy
property
higher_order_energy
Fits a quadratic (order=2) to the Energy vs. time, and returns the energy with this subtracted.
higher_order_energy_spread
property
higher_order_energy_spread
Legacy syntax to compute the standard deviation of higher_order_energy.
id
property
writable
id
id integer
in_t_coordinates
property
in_t_coordinates
Returns True if all particles have the same t coordinate
in_z_coordinates
property
in_z_coordinates
Returns True if all particles have the same z coordinate
kinetic_energy
property
kinetic_energy
Kinetic energy in eV
mass
property
mass
Rest mass in eV
n_alive
property
n_alive
Number of alive particles, defined by status == 1
n_dead
property
n_dead
Number of alive particles, defined by status != 1
n_particle
property
n_particle
Total number of particles. Same as len
norm_emit_4d
property
norm_emit_4d
Normalized emittance in the xy planes (4D)
norm_emit_x
property
norm_emit_x
Normalized emittance in the x plane
norm_emit_y
property
norm_emit_y
Normalized emittance in the x plane
p
property
p
Total momemtum in eV/c
pr
property
pr
Momentum in the radial direction in eV/c r_hat = cos(theta) xhat + sin(theta) yhat pr = p dot r_hat
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
px
property
writable
px
px coordinate in [eV/c]
px_bar
property
px_bar
Normalized px in units of sqrt(m)
py
property
writable
py
py coordinate in [eV/c]
py_bar
property
py_bar
Normalized py in units of sqrt(m)
pz
property
writable
pz
pz coordinate in [eV/c]
r
property
r
Radius in the xy plane: r = sqrt(x^2 + y^2) in m
species
property
writable
species
species string
species_charge
property
species_charge
Species charge in C
status
property
writable
status
status coordinate in [1]
t
property
writable
t
t coordinate in [s]
theta
property
theta
Angle in xy plane: theta = arctan2(y, x) in radians
weight
property
writable
weight
weight coordinate in [C]
x
property
writable
x
x coordinate in [m]
x_bar
property
x_bar
Normalized x in units of sqrt(m)
xp
property
xp
x slope px/pz (dimensionless)
y
property
writable
y
y coordinate in [m]
y_bar
property
y_bar
Normalized y in units of sqrt(m)
yp
property
yp
y slope py/pz (dimensionless)
z
property
writable
z
z coordinate in [m]
__add__
__add__(other)
Overloads the + operator to join particle groups. Simply calls join_particle_groups
Source code in beamphysics/particles.py
1823 1824 1825 1826 1827 1828 | |
__contains__
__contains__(item)
Checks internal data
Source code in beamphysics/particles.py
1830 1831 1832 | |
__eq__
__eq__(other)
Check equality of internal data
Source code in beamphysics/particles.py
1834 1835 1836 1837 1838 1839 1840 1841 1842 | |
__getitem__
__getitem__(key: str)
Returns a property or statistical quantity that can be computed:
P['x']returns the x arrayP['sigmx_x']returns the std(x) scalarP['norm_emit_x']returns the norm_emit_x scalar
Parts can also be given. Example: P[0:10] returns a new ParticleGroup with the first 10 elements.
Source code in beamphysics/particles.py
894 895 896 897 898 899 900 901 902 903 904 905 906 907 908 909 910 911 912 913 914 915 916 917 918 919 920 921 922 923 924 925 926 927 928 929 930 931 932 933 934 935 936 937 938 939 940 941 942 943 944 945 946 947 948 949 950 951 952 953 954 955 | |
__len__
__len__()
Return the number of particles.
Source code in beamphysics/particles.py
1844 1845 1846 | |
__repr__
__repr__()
Return a detailed string representation with memory location.
Source code in beamphysics/particles.py
1853 1854 1855 1856 | |
__str__
__str__()
Return a human-readable summary string.
Source code in beamphysics/particles.py
1848 1849 1850 1851 | |
apply_wakefield
apply_wakefield(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, default:False) –If True, modifies in place. If False, returns a modified copy. Default is False.
-
include_self_kick(bool, default:True) –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)
Source code in beamphysics/particles.py
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 1597 1598 1599 1600 1601 1602 1603 1604 1605 1606 1607 1608 1609 1610 1611 1612 1613 1614 1615 | |
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
382 383 384 385 386 387 388 389 390 | |
avg
avg(key: str)
Return the weighted average of a particle property.
Parameters:
-
key(str) –Any valid particle property name.
Returns:
-
float–Weighted average using
self.weightas weights.
Source code in beamphysics/particles.py
689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 | |
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:
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
wavelengthis not a positive number.
Source code in beamphysics/particles.py
833 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 | |
copy
copy()
Return a deep copy of this ParticleGroup.
Returns:
-
ParticleGroup–An independent copy with no shared data.
Source code in beamphysics/particles.py
1766 1767 1768 1769 1770 1771 1772 1773 1774 1775 | |
cov
cov(*keys)
Covariance matrix from any properties
Example: P = ParticleGroup(h5) P.cov('x', 'px', 'y', 'py')
Source code in beamphysics/particles.py
728 729 730 731 732 733 734 735 736 737 738 | |
delta
delta(key: str)
Return an attribute array relative to its weighted mean.
Parameters:
-
key(str) –Any valid particle property name.
Returns:
-
ndarray–Array of
self[key] - self.avg(key).
Source code in beamphysics/particles.py
623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 | |
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:
-
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.
Source code in beamphysics/particles.py
985 986 987 988 989 990 991 992 993 994 995 996 997 998 999 1000 | |
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
1023 1024 1025 1026 1027 1028 1029 1030 1031 1032 1033 1034 | |
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:
-
z(float, default:None) –Target z position in meters. If None, particles are drifted to the mean z.
Source code in beamphysics/particles.py
1002 1003 1004 1005 1006 1007 1008 1009 1010 1011 1012 1013 1014 1015 1016 1017 1018 1019 1020 1021 | |
fractional_split
fractional_split(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.
Source code in beamphysics/particles.py
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 1746 1747 1748 1749 1750 1751 1752 1753 1754 1755 1756 1757 1758 1759 1760 1761 1762 1763 1764 | |
from_bmad
classmethod
from_bmad(bmad_dict)
Convert Bmad particle data to a ParticleGroup.
This reverses the conversion done by to_bmad, mapping Bmad phase space coordinates back to a ParticleGroup data format.
Parameters:
-
bmad_dict(dict) –A dictionary containing Bmad phase space coordinates with keys: 'x', 'px', 'y', 'py', 'z', 'pz', 'charge', 'species', 'tref', 'state'
Returns:
Source code in beamphysics/particles.py
1073 1074 1075 1076 1077 1078 1079 1080 1081 1082 1083 1084 1085 1086 1087 1088 1089 1090 1091 1092 1093 1094 | |
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
484 485 486 487 488 489 490 491 492 493 494 495 496 497 498 499 | |
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
740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 | |
info
info(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
Source code in beamphysics/particles.py
875 876 877 878 879 880 881 882 883 884 885 886 887 888 889 890 891 892 | |
linear_point_transform
linear_point_transform(mat3: Union[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
mat3is 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
1860 1861 1862 1863 1864 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 | |
max
max(key: str)
Return the maximum value of a particle property.
Parameters:
-
key(str) –Any valid particle property name.
Returns:
-
float–Maximum value.
Source code in beamphysics/particles.py
657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 | |
min
min(key: str)
Return the minimum value of a particle property.
Parameters:
-
key(str) –Any valid particle property name.
Returns:
-
float–Minimum value.
Source code in beamphysics/particles.py
641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 | |
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:
-
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
Source code in beamphysics/particles.py
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 1440 1441 1442 1443 1444 1445 1446 1447 1448 1449 1450 1451 1452 1453 1454 1455 1456 1457 1458 | |
ptp
ptp(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.
Source code in beamphysics/particles.py
673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 | |
resample
resample(n=0, equal_weights=False)
Resample particles randomly.
If n equals self.n_particle or n=0, particle indices will be scrambled.
Otherwise if weights are equal, a random subset of particles will be selected.
Otherwise if weights are not equal, particles will be sampled according to their weight using scipy.stats.rv_discrete. Note that this latter method can result in duplicate particles, and can be very slow for a large number of particles.
Parameters:
-
n(int, default:0) –Number to resample. If n=0, this will use all particles.
-
equal_weights(bool, default:False) –If True, will ensure that all particles have equal weights.
Returns:
Source code in beamphysics/particles.py
1777 1778 1779 1780 1781 1782 1783 1784 1785 1786 1787 1788 1789 1790 1791 1792 1793 1794 1795 1796 1797 1798 1799 1800 1801 1802 1803 1804 1805 1806 | |
rotate
rotate(*, 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, default:0.0) –Rotation around the x axis, radians
-
y_rot(float, default:0.0) –Rotation around the y axis, radians
-
z_rot(float, default:0.0) –Rotation around the z axis, radians
-
order(str, default:'zxy') –A 3-character string specifying the rotation order (e.g., 'zxy', 'zyx'). Each character must be one of 'x', 'y', or 'z'.
-
xc(float, default:0.0) –Center of rotation, x coordinate. Default to zero
-
yc(float, default:0.0) –Center of rotation, x coordinate. Default to zero
-
zc(float, default:0.0) –Center of rotation, x coordinate. Default to zero
Notes
This method modifies the object in-place.
Source code in beamphysics/particles.py
1910 1911 1912 1913 1914 1915 1916 1917 1918 1919 1920 1921 1922 1923 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 1962 1963 1964 1965 1966 1967 | |
rotate_x
rotate_x(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, default:0.0) –Center of rotation, y coordinate. Default to zero
-
zc(float, default:0.0) –Center of rotation, z coordinate. Default to zero
Source code in beamphysics/particles.py
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 2001 2002 2003 2004 2005 2006 | |
rotate_y
rotate_y(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, default:0.0) –Center of rotation, x coordinate. Default to zero
-
zc(float, default:0.0) –Center of rotation, z coordinate. Default to zero
Source code in beamphysics/particles.py
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 2040 2041 2042 2043 2044 2045 | |
rotate_z
rotate_z(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, default:0.0) –Center of rotation, x coordinate. Default to zero
-
yc(float, default:0.0) –Center of rotation, y coordinate. Default to zero
Source code in beamphysics/particles.py
2047 2048 2049 2050 2051 2052 2053 2054 2055 2056 2057 2058 2059 2060 2061 2062 2063 2064 2065 2066 2067 2068 2069 2070 2071 2072 2073 2074 2075 2076 2077 2078 2079 2080 2081 2082 2083 2084 | |
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:
-
*keys(str, default:()) –Keys to plot slice statistics for (e.g.,
'sigma_x','norm_emit_x'). -
n_slice(int, default:100) –Number of slices. Default is 100.
-
slice_key(str, default:None) –Key to slice on. If None, defaults to
'z'for t-coordinates or't'for z-coordinates. -
tex(bool, default:True) –Use TeX for axis labels. Default is True.
-
nice(bool, default:True) –Scale to nice units. Default is True.
-
return_figure(bool, default:False) –If True, return the matplotlib Figure. Default is False.
-
xlim(tuple of float, default:None) –Manual x-axis limits in raw units.
-
ylim(tuple of float, default:None) –Manual y-axis limits in raw units.
-
**kwargs–Additional keyword arguments passed to
plt.subplots.
Returns:
-
None or Figure–Returns a Figure only if
return_figure=True.
Source code in beamphysics/particles.py
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 1544 1545 1546 1547 1548 1549 1550 1551 1552 1553 1554 1555 1556 1557 1558 1559 1560 1561 1562 | |
slice_statistics
slice_statistics(*keys, n_slice=100, slice_key=None)
Compute slice statistics along the bunch.
Parameters:
-
*keys(str, default:()) –Statistical quantities to compute per slice (e.g.,
'sigma_x','norm_emit_x','mean_energy'). -
n_slice(int, default:100) –Number of slices. Default is 100.
-
slice_key(str, default:None) –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.
Source code in beamphysics/particles.py
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 1487 1488 1489 1490 1491 1492 1493 1494 1495 1496 1497 1498 1499 1500 1501 1502 1503 1504 1505 | |
split
split(n_chunks=100, key='z')
Split the particle group into even chunks sorted by a key.
Parameters:
-
n_chunks(int, default:100) –Number of chunks to split into. Default is 100.
-
key(str, default:'z') –Particle attribute to sort by before splitting. Default is
'z'.
Returns:
-
list of ParticleGroup–A list of ParticleGroup objects, one per chunk.
Source code in beamphysics/particles.py
1688 1689 1690 1691 1692 1693 1694 1695 1696 1697 1698 1699 1700 1701 1702 1703 1704 | |
std
std(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.weightas weights.
Source code in beamphysics/particles.py
708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 | |
to_bmad
to_bmad(p0c=None, tref=None)
Convert to Bmad phase space coordinates.
Bmad openPMD-beamphysics
Bmad x = x Bmad px = px/p0c Bmad y = y Bmad py = py/p0c Bmad z = -beta * c * (t - tref) Bmad pz = p/p0c - 1 Bmad t = t
Parameters:
-
p0c(float, default:None) –Reference momentum times the speed of light (in eV). Default is None, which uses self['mean_p'].
-
tref(float, default:None) –Reference time (in seconds). Default is None, which uses self['mean_t'].
Returns:
-
dict–A dictionary containing Bmad phase space coordinates.
Source code in beamphysics/particles.py
1043 1044 1045 1046 1047 1048 1049 1050 1051 1052 1053 1054 1055 1056 1057 1058 1059 1060 1061 1062 1063 1064 1065 1066 1067 1068 1069 1070 1071 | |
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
777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 | |
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
794 795 796 797 798 799 800 801 802 803 | |
units
units(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–
Source code in beamphysics/particles.py
407 408 409 410 411 412 413 414 415 416 417 418 419 420 | |
wakefield_plot
wakefield_plot(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, default:None) –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(Axes, default:None) –An existing Axes to plot into. If None, a new figure and axes are created.
-
xlim(tuple of float, default:None) –Limits to apply to the x-axis, in native units.
-
ylim(tuple of float, default:None) –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, default:None) –Number of bins to use for the density histogram.
-
kwargs(dict, default:{}) –Additional keyword arguments passed to
plt.subplots()if a new axis is created.
Returns:
-
fig(Figure) –The matplotlib figure containing the plot.
Source code in beamphysics/particles.py
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 1667 1668 1669 1670 1671 1672 1673 1674 1675 1676 1677 1678 1679 1680 1681 1682 1683 1684 1685 | |
where
where(x)
Filter particles by a boolean condition.
Parameters:
-
x(array_like of bool) –Boolean mask or condition array. Particles where
xis True are selected.
Returns:
-
ParticleGroup–A new ParticleGroup containing only the selected particles.
Source code in beamphysics/particles.py
957 958 959 960 961 962 963 964 965 966 967 968 969 970 971 972 | |
write
write(h5, name=None) -> None
Write particle data to an HDF5 file or group in openPMD format.
Parameters:
-
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. -
name(str, default:None) –Name for the subgroup/bunch written inside the "particles" group (or provided group). If None,
write_pmd_bunchwill 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")
Source code in beamphysics/particles.py
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 1347 1348 1349 1350 1351 1352 1353 1354 1355 1356 1357 1358 1359 1360 1361 1362 1363 1364 1365 1366 1367 1368 | |
write_astra
write_astra(filePath, verbose=False, probe=False)
Write Astra style particles.
For now, the species must be electrons.
If probe, the six standard probe particles will be written.
Source code in beamphysics/particles.py
1099 1100 1101 1102 1103 1104 1105 1106 1107 | |
write_bmad
write_bmad(filePath, p0c=None, t_ref=0, verbose=False)
Write particle data in Bmad format.
Parameters:
-
filePath(str or Path) –Output file path.
-
p0c(float, default:None) –Reference momentum in eV/c. If None, the average pz is used.
-
t_ref(float, default:0) –Reference time in seconds. Default is 0.
-
verbose(bool, default:False) –If True, print information during writing. Default is False.
Source code in beamphysics/particles.py
1109 1110 1111 1112 1113 1114 1115 1116 1117 1118 1119 1120 1121 1122 1123 1124 | |
write_elegant
write_elegant(filePath, verbose=False)
Write particle data in Elegant SDDS format.
Parameters:
-
filePath(str or Path) –Output file path.
-
verbose(bool, default:False) –If True, print information during writing. Default is False.
Source code in beamphysics/particles.py
1126 1127 1128 1129 1130 1131 1132 1133 1134 1135 1136 1137 | |
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:
-
filePath(str or Path) –Output file path.
-
n_slice(int, default:None) –Number of slices. If None, a default is determined automatically.
-
verbose(bool, default:False) –If True, print information during writing. Default is False.
Source code in beamphysics/particles.py
1139 1140 1141 1142 1143 1144 1145 1146 1147 1148 1149 1150 1151 1152 1153 1154 1155 | |
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:
-
filePath(str or Path) –Output file path.
-
n_slice(int, default:None) –Number of slices. If None, a default is determined automatically.
-
return_input_str(bool, default:False) –Defaults to False.
-
verbose(bool, default:False) –If True, print information during writing. Default is False.
Source code in beamphysics/particles.py
1157 1158 1159 1160 1161 1162 1163 1164 1165 1166 1167 1168 1169 1170 1171 1172 1173 1174 1175 1176 1177 1178 1179 1180 | |
write_genesis4_distribution
write_genesis4_distribution(filePath, verbose=False)
Write particle data as a Genesis 1.3 v4 distribution HDF5 file.
Parameters:
-
filePath(str or Path) –Output file path.
-
verbose(bool, default:False) –If True, print information during writing. Default is False.
Source code in beamphysics/particles.py
1182 1183 1184 1185 1186 1187 1188 1189 1190 1191 1192 1193 | |
write_gpt
write_gpt(filePath, asci2gdf_bin=None, verbose=False)
Write particle data in GPT (General Particle Tracer) format.
Parameters:
-
filePath(str or Path) –Output file path.
-
asci2gdf_bin(str, default:None) –Path to the
asci2gdfbinary for GDF conversion. If None, an ASCII file is written. -
verbose(bool, default:False) –If True, print information during writing. Default is False.
Source code in beamphysics/particles.py
1195 1196 1197 1198 1199 1200 1201 1202 1203 1204 1205 1206 1207 1208 1209 | |
write_impact
write_impact(filePath, cathode_kinetic_energy_ref=None, include_header=True, verbose=False)
Write particle data in Impact-T format.
Parameters:
-
filePath(str or Path) –Output file path.
-
cathode_kinetic_energy_ref(float, default:None) –Reference kinetic energy at the cathode in eV. If None, the average kinetic energy is used.
-
include_header(bool, default:True) –If True, include the Impact-T header line. Default is True.
-
verbose(bool, default:False) –If True, print information during writing. Default is False.
Source code in beamphysics/particles.py
1211 1212 1213 1214 1215 1216 1217 1218 1219 1220 1221 1222 1223 1224 1225 1226 1227 1228 1229 1230 1231 1232 1233 1234 1235 1236 1237 1238 1239 1240 1241 | |
write_litrack
write_litrack(filePath, p0c=None, verbose=False)
Write particle data in LiTrack format.
Parameters:
-
filePath(str or Path) –Output file path.
-
p0c(float, default:None) –Reference momentum in eV/c. If None, the average pz is used.
-
verbose(bool, default:False) –If True, print information during writing. Default is False.
Source code in beamphysics/particles.py
1243 1244 1245 1246 1247 1248 1249 1250 1251 1252 1253 1254 1255 1256 | |
write_lucretia
write_lucretia(filePath, ele_name='BEGINNING', t_ref=0, stop_ix=None, verbose=False)
Write particle data in Lucretia format.
Parameters:
-
filePath(str or Path) –Output file path.
-
ele_name(str, default:'BEGINNING') –Element name for the Lucretia beam. Default is
'BEGINNING'. -
t_ref(float, default:0) –Reference time in seconds. Default is 0.
-
stop_ix(int, default:None) –Stop index for the Lucretia beam structure.
-
verbose(bool, default:False) –If True, print information during writing. Default is False.
Source code in beamphysics/particles.py
1258 1259 1260 1261 1262 1263 1264 1265 1266 1267 1268 1269 1270 1271 1272 1273 1274 1275 1276 1277 1278 1279 | |
write_opal
write_opal(filePath, verbose=False, dist_type='emitted')
Write particle data in OPAL (Object Oriented Parallel Accelerator Library) format.
Parameters:
-
filePath(str or Path) –Output file path.
-
verbose(bool, default:False) –If True, print information during writing. Default is False.
-
dist_type(str, default:'emitted') –Distribution type, either
'emitted'or'injected'. Default is'emitted'.
Source code in beamphysics/particles.py
1301 1302 1303 1304 1305 1306 1307 1308 1309 1310 1311 1312 1313 1314 1315 | |
write_simion
write_simion(filePath, color=0, flip_z_to_x=True, verbose=False)
Write particle data in SIMION format.
Parameters:
-
filePath(str or Path) –Output file path.
-
color(int, default:0) –SIMION color index for particles. Default is 0.
-
flip_z_to_x(bool, default:True) –If True, swap z and x axes for SIMION's coordinate convention. Default is True.
-
verbose(bool, default:False) –If True, print information during writing. Default is False.
Source code in beamphysics/particles.py
1281 1282 1283 1284 1285 1286 1287 1288 1289 1290 1291 1292 1293 1294 1295 1296 1297 1298 1299 | |