Xfields API

Beam elements

Space Charge 3D

class xfields.SpaceCharge3D(_context=None, _buffer=None, _offset=None, update_on_track=True, length=None, apply_z_kick=True, fieldmap=None, x_range=None, y_range=None, z_range=None, nx=None, ny=None, nz=None, dx=None, dy=None, dz=None, x_grid=None, y_grid=None, z_grid=None, rho=None, phi=None, solver=None, gamma0=None, fftplan=None)

Simulates the effect of space charge on a bunch.

Parameters
  • context (XfContext) – identifies the context on which the computation is executed.

  • update_on_track (bool) – If True the beam field map is update at each interaction. If False the initial field map is used at each interaction (frozen model). The default is True.

  • length (float) – the length of the space-charge interaction in meters.

  • apply_z_kick (bool) – If True, the longitudinal kick on the particles is applied.

  • x_range (tuple) – Horizontal extent (in meters) of the computing grid.

  • y_range (tuple) – Vertical extent (in meters) of the computing grid.

  • z_range (tuple) – Longitudina extent (in meters) of the computing grid.

  • nx (int) – Number of cells in the horizontal direction.

  • ny (int) – Number of cells in the vertical direction.

  • nz (int) – Number of cells in the vertical direction.

  • dx (float) – Horizontal cell size in meters. It can be provided alternatively to nx.

  • dy (float) – Vertical cell size in meters. It can be provided alternatively to ny.

  • dz (float) – Longitudinal cell size in meters.It can be provided alternatively to nz.

  • x_grid (np.ndarray) – Equispaced array with the horizontal grid points (cell centers). It can be provided alternatively to x_range, dx/nx.

  • y_grid (np.ndarray) – Equispaced array with the horizontal grid points (cell centers). It can be provided alternatively to y_range, dy/ny.

  • z_grid (np.ndarray) – Equispaced array with the horizontal grid points (cell centers). It can be provided alternatively to z_range, dz/nz.

  • rho (np.ndarray) – initial charge density at the grid points in Coulomb/m^3.

  • phi (np.ndarray) – initial electric potential at the grid points in Volts. If not provided the phi is calculated from rho using the Poisson solver (if available).

  • solver (str or solver object) – Defines the Poisson solver to be used to compute phi from rho. Accepted values are FFTSolver3D and FFTSolver2p5D. A Xfields solver object can also be provided. In case update_on_track``is ``False and phi is provided by the user, this argument can be omitted.

  • gamma0 (float) – Relativistic gamma factor of the beam. This is required only if the solver is FFTSolver3D.

Returns

A space-charge 3D beam element.

Return type

(SpaceCharge3D)

copy(_context=None, _buffer=None, _offset=None)
property iscollective
track(particles)

Computes and applies the space-charge forces for the provided set of particles.

Parameters

particles (Particles Object) – Particles to be tracked.

fieldmap
length

Space Charge Bi-Gaussian

class xfields.SpaceChargeBiGaussian(_context=None, _buffer=None, _offset=None, _xobject=None, update_on_track=False, length=None, apply_z_kick=False, longitudinal_profile=None, mean_x=0.0, mean_y=0.0, sigma_x=None, sigma_y=None, fieldmap=None, min_sigma_diff=1e-10, **kwargs)
to_dict()
length
longitudinal_profile
fieldmap
track(particles)
property iscollective
property mean_x
property mean_y
property sigma_x
property sigma_y

Beam-beam Bi-Gaussian 2D

class xfields.BeamBeamBiGaussian2D(scale_strength=1.0, other_beam_q0=None, other_beam_beta0=None, other_beam_num_particles=None, other_beam_Sigma_11=None, other_beam_Sigma_13=None, other_beam_Sigma_33=None, ref_shift_x=0, ref_shift_y=0, other_beam_shift_x=0, other_beam_shift_y=0, post_subtract_px=0, post_subtract_py=0, min_sigma_diff=1e-10, config_for_update=None, **kwargs)
other_beam_num_particles
other_beam_Sigma_11
other_beam_Sigma_13
other_beam_Sigma_33
other_beam_q0
other_beam_beta0
ref_shift_x
ref_shift_y
other_beam_shift_x
other_beam_shift_y
post_subtract_px
post_subtract_py
min_sigma_diff
scale_strength
property n_particles
property q0
property beta0
property mean_x
property mean_y
property sigma_x
property sigma_y
property d_px
property d_py

Beam-beam Bi-Gaussian 3D

class xfields.BeamBeamBiGaussian3D(phi=None, alpha=None, other_beam_q0=None, scale_strength=1.0, slices_other_beam_num_particles=None, slices_other_beam_x_center=0.0, slices_other_beam_px_center=0.0, slices_other_beam_y_center=0.0, slices_other_beam_py_center=0.0, slices_other_beam_zeta_center=None, slices_other_beam_pzeta_center=0.0, flag_beamstrahlung=0, slices_other_beam_zeta_bin_width_star_beamstrahlung=None, other_beam_sigma_55_star_beamstrahlung=None, slices_other_beam_x_center_star=None, slices_other_beam_px_center_star=None, slices_other_beam_y_center_star=None, slices_other_beam_py_center_star=None, slices_other_beam_zeta_center_star=None, slices_other_beam_pzeta_center_star=None, slices_other_beam_Sigma_11=None, slices_other_beam_Sigma_12=None, slices_other_beam_Sigma_13=None, slices_other_beam_Sigma_14=None, slices_other_beam_Sigma_22=None, slices_other_beam_Sigma_23=None, slices_other_beam_Sigma_24=None, slices_other_beam_Sigma_33=None, slices_other_beam_Sigma_34=None, slices_other_beam_Sigma_44=None, slices_other_beam_Sigma_11_star=None, slices_other_beam_Sigma_12_star=None, slices_other_beam_Sigma_13_star=None, slices_other_beam_Sigma_14_star=None, slices_other_beam_Sigma_22_star=None, slices_other_beam_Sigma_23_star=None, slices_other_beam_Sigma_24_star=None, slices_other_beam_Sigma_33_star=None, slices_other_beam_Sigma_34_star=None, slices_other_beam_Sigma_44_star=None, ref_shift_x=0, ref_shift_px=0, ref_shift_y=0, ref_shift_py=0, ref_shift_zeta=0, ref_shift_pzeta=0, other_beam_shift_x=0, other_beam_shift_px=0, other_beam_shift_y=0, other_beam_shift_py=0, other_beam_shift_zeta=0, other_beam_shift_pzeta=0, post_subtract_x=0, post_subtract_px=0, post_subtract_y=0, post_subtract_py=0, post_subtract_zeta=0, post_subtract_pzeta=0, min_sigma_diff=1e-10, threshold_singular=1e-28, old_interface=None, config_for_update=None, _sin_phi=None, _cos_phi=None, _tan_phi=None, _sin_alpha=None, _cos_alpha=None, **kwargs)
num_slices_other_beam
slices_other_beam_num_particles
other_beam_q0
scale_strength
ref_shift_x
ref_shift_px
ref_shift_y
ref_shift_py
ref_shift_zeta
ref_shift_pzeta
other_beam_shift_x
other_beam_shift_px
other_beam_shift_y
other_beam_shift_py
other_beam_shift_zeta
other_beam_shift_pzeta
post_subtract_x
post_subtract_px
post_subtract_y
post_subtract_py
post_subtract_zeta
post_subtract_pzeta
min_sigma_diff
threshold_singular
slices_other_beam_zeta_bin_width_star_beamstrahlung
other_beam_sigma_55_star_beamstrahlung
property flag_beamstrahlung
update_from_recieved_moments()
property sin_phi
property cos_phi
property tan_phi
property sin_alpha
property cos_alpha
property phi
property alpha
property slices_other_beam_x_center
property slices_other_beam_px_center
property slices_other_beam_y_center
property slices_other_beam_py_center
property slices_other_beam_zeta_center
property slices_other_beam_pzeta_center
property slices_other_beam_Sigma_11
property slices_other_beam_Sigma_12
property slices_other_beam_Sigma_13
property slices_other_beam_Sigma_14
property slices_other_beam_Sigma_22
change_back_ref_frame_and_subtract_dipolar
change_ref_frame
slices_other_beam_Sigma_11_star
slices_other_beam_Sigma_12_star
slices_other_beam_Sigma_13_star
slices_other_beam_Sigma_14_star
slices_other_beam_Sigma_22_star
property slices_other_beam_Sigma_23
slices_other_beam_Sigma_23_star
slices_other_beam_Sigma_24_star
slices_other_beam_Sigma_33_star
slices_other_beam_Sigma_34_star
slices_other_beam_Sigma_44_star
slices_other_beam_px_center_star
slices_other_beam_py_center_star
slices_other_beam_pzeta_center_star
slices_other_beam_x_center_star
slices_other_beam_y_center_star
slices_other_beam_zeta_center_star
synchro_beam_kick
property slices_other_beam_Sigma_24
property slices_other_beam_Sigma_33
property slices_other_beam_Sigma_34
property slices_other_beam_Sigma_44

Field Maps

Beam-beam Tri-linear Iterpolated Field Map

class xfields.fieldmaps.TriLinearInterpolatedFieldMap(_context=None, _buffer=None, _offset=None, _xobject=None, x_range=None, y_range=None, z_range=None, nx=None, ny=None, nz=None, dx=None, dy=None, dz=None, x_grid=None, y_grid=None, z_grid=None, rho=None, phi=None, solver=None, scale_coordinates_in_solver=(1.0, 1.0, 1.0), updatable=True, fftplan=None)

Builds a linear interpolator for a 3D field map. The map can be updated using the Parcle In Cell method.

Parameters
  • context (xobjects context) – identifies the context on which the computation is executed.

  • x_range (tuple) – Horizontal extent (in meters) of the computing grid.

  • y_range (tuple) – Vertical extent (in meters) of the computing grid.

  • z_range (tuple) – Longitudina extent (in meters) of the computing grid.

  • nx (int) – Number of cells in the horizontal direction.

  • ny (int) – Number of cells in the vertical direction.

  • nz (int) – Number of cells in the vertical direction.

  • dx (float) – Horizontal cell size in meters. It can be provided alternatively to nx.

  • dy (float) – Vertical cell size in meters. It can be provided alternatively to ny.

  • dz (float) – Longitudinal cell size in meters.It can be provided alternatively to nz.

  • x_grid (np.ndarray) – Equispaced array with the horizontal grid points (cell centers). It can be provided alternatively to x_range, dx/nx.

  • y_grid (np.ndarray) – Equispaced array with the horizontal grid points (cell centers). It can be provided alternatively to y_range, dy/ny.

  • z_grid (np.ndarray) – Equispaced array with the horizontal grid points (cell centers). It can be provided alternatively to z_range, dz/nz.

  • rho (np.ndarray) – initial charge density at the grid points in Coulomb/m^3.

  • phi (np.ndarray) – initial electric potential at the grid points in Volts. If not provided the phi is calculated from rho using the Poisson solver (if available).

  • solver (str or solver object) – Defines the Poisson solver to be used to compute phi from rho. Accepted values are FFTSolver3D and FFTSolver2p5D. A Xfields solver object can also be provided. In case update_on_track``is ``False and phi is provided by the user, this argument can be omitted.

  • scale_coordinates_in_solver (tuple) – Three coefficients used to rescale the grid coordinates in the definition of the solver. The default is (1.,1.,1.).

  • updatable (bool) – If True the field map can be updated after creation. Default is True.

Returns

Interpolator object.

Return type

(TriLinearInterpolatedFieldMap)

get_values_at_points(x, y, z, return_rho=True, return_phi=True, return_dphi_dx=True, return_dphi_dy=True, return_dphi_dz=True)

Returns the charge density, the field potential and its derivatives at the points specified by x, y, z. The output can be customized (see below). Zeros are returned for points outside the grid.

Parameters
  • x (float64 array) – Horizontal coordinates at which the field is evaluated.

  • y (float64 array) – Vertical coordinates at which the field is evaluated.

  • z (float64 array) – Longitudinal coordinates at which the field is evaluated.

  • return_rho (bool) – If True, the charge density at the given points is returned.

  • return_phi (bool) – If True, the potential at the given points is returned.

  • return_dphi_dx (bool) – If True, the horizontal derivative of the potential at the given points is returned.

  • return_dphi_dy – If True, the vertical derivative of the potential at the given points is returned.

  • return_dphi_dz – If True, the longitudinal derivative of the potential at the given points is returned.

Returns

The required quantities at the provided points.

Return type

(tuple of float64 array)

update_from_particles(particles=None, x_p=None, y_p=None, z_p=None, ncharges_p=None, state_p=None, q0_coulomb=None, reset=True, update_phi=True, solver=None, force=False)

Updates the charge density at the grid using a given set of particles, which can be provided by a particles object or by individual arrays. The potential can be optionally updated accordingly.

Parameters
  • particles (xtrack.Particles) – xtrack particle object.

  • x_p (float64 array) – Horizontal coordinates of the macroparticles.

  • y_p (float64 array) – Vertical coordinates of the macroparticles.

  • z_p (float64 array) – Longitudinal coordinates of the macroparticles.

  • ncharges_p (float64 array) – Number of reference charges in the macroparticles.

  • state_p (int64, array) – particle state (>0 active, lost otherwise)

  • q0_coulomb (float64) – Reference charge in Coulomb.

  • reset (bool) – If True the stored charge density is overwritten with the provided one. If False the provided charge density is added to the stored one. The default is True.

  • update_phi (bool) – If True the stored potential is recalculated from the stored charge density.

  • solver (Solver object) – solver object to be used to solve Poisson’s equation (compute phi from rho). If None is provided the solver attached to the fieldmap is used (if any). The default is None.

  • force (bool) – If True the potential is updated even if the map is declared as not updateable. The default is False.

update_rho(rho, reset=True, force=False)

Updates the charge density on the grid.

Parameters
  • rho (float64 array) – Charge density at the grid points in C/m^3.

  • reset (bool) – If True the stored charge density is overwritten with the provided one. If False the provided charge density is added to the stored one. The default is True.

  • force (bool) – If True the charge density is updated even if the map is declared as not updateable. The default is False.

update_phi(phi, reset=True, force=False)

Updates the potential on the grid. The stored derivatives are also updated.

Parameters
  • rho (float64 array) – Potential at the grid points.

  • reset (bool) – If True the stored potential is overwritten with the provided one. If False the provided potential is added to the stored one. The default is True.

  • force (bool) – If True the potential is updated even if the map is declared as not updateable. The default is False.

update_phi_from_rho(solver=None)

Updates the potential on the grid (phi) from the charge density on the grid (phi). It requires a Poisson solver object. If none is provided the one attached to the fieldmap is used (if any).

Parameters

solver (Solver object) – solver object to be used to solve Poisson’s equation. If None is provided the solver attached to the fieldmap is used (if any). The default is None.

generate_solver(solver, fftplan)

Generates a Poisson solver associated to the defined grid.

Parameters
  • solver (str) – Defines the Poisson solver to be used

  • and (to compute phi from rho. Accepted values are FFTSolver3D) –

  • FFTSolver2p5D.

Returns

Solver object associated to the defined grid.

Return type

(Solver)

property x_grid

Array with the horizontal grid points (cell centers).

property y_grid

Array with the vertical grid points (cell centers).

property z_grid

Array with the longitudinal grid points (cell centers).

property nx

Number of cells in the horizontal direction.

property ny

Number of cells in the vertical direction.

property nz

Number of cells in the longitudinal direction.

property dx

Horizontal cell size in meters.

property dy

Vertical cell size in meters.

property dz

Longitudinal cell size in meters.

property rho
property phi

Electric potential at the grid points in Volts.

property dphi_dx
property dphi_dy
property dphi_dz

Beam-beam Bi-Gaussian Field Map

class xfields.fieldmaps.BiGaussianFieldMap(_context=None, _buffer=None, _offset=None, _xobject=None, mean_x=0.0, mean_y=0.0, sigma_x=None, sigma_y=None, min_sigma_diff=1e-10, updatable=True)

Builds a field-map object that provides the fields generated by a normalized 2D Gaussian distribution (Bassetti-Erskine formula).

Parameters
  • context (xobjects context) – identifies the context on which the computation is executed.

  • mean_x (float64) – Horizontal position (in meters) of the Gaussian distribution. It can be updated after the object creation. Default is 0..

  • mean_y (float64) – Vertical position (in meters) of the Gaussian distribution. It can be updated after the object creation. Default is 0..

  • sigma_x (float64) – Horizontal r.m.s. size (in meters) of the Gaussian distribution. It can be updated after the object creation. Default is None.

  • sigma_y (float64) – Vertical r.m.s. size (in meters) of tthe Gaussian distribution. It can be updated after the object creation. Default is None.

  • min_sigma_diff (float64) – Difference between sigma_x and sigma_y (in meters) below which round distribution is assumed.

  • updatable (bool) – If True the field map can be updated after creation. Default is True.

Returns

Field map object.

Return type

(BiGaussianFieldMap)

mean_x
mean_y
sigma_x
sigma_y
min_sigma_diff
property updatable
update_from_particles(x_p, y_p, z_p, ncharges_p, q0_coulomb, reset=True, update_phi=True, solver=None, force=False)
update_rho(rho, reset)
update_phi(phi, reset=True, force=False)
update_phi_from_rho(solver=None)
generate_solver(solver)

Solvers

FFT Solver 3D

class xfields.solvers.FFTSolver3D(dx, dy, dz, nx, ny, nz, context=None, fftplan=None)

Creates a Poisson solver object that solves the full 3D Poisson equation using the FFT method (free space).

Parameters
  • nx (int) – Number of cells in the horizontal direction.

  • ny (int) – Number of cells in the vertical direction.

  • nz (int) – Number of cells in the vertical direction.

  • dx (float) – Horizontal cell size in meters.

  • dy (float) – Vertical cell size in meters.

  • dz (float) – Longitudinal cell size in meters.

  • context (XfContext) – identifies the context on which the computation is executed.

Returns

Poisson solver object.

Return type

(FFTSolver3D)

solve(rho)

Solves Poisson’s equation in free space for a given charge density.

Parameters

rho (float64 array) – charge density at the grid points in Coulomb/m^3.

Returns

electric potential at the grid points in Volts.

Return type

phi (float64 array)

FFT Solver 2.5D

class xfields.solvers.FFTSolver2p5D(dx, dy, dz, nx, ny, nz, context=None, fftplan=None)

Creates a Poisson solver object that solve’s Poisson equation in the 2.5D aaoroximation equation the FFT method (free space).

Parameters
  • nx (int) – Number of cells in the horizontal direction.

  • ny (int) – Number of cells in the vertical direction.

  • nz (int) – Number of cells in the vertical direction.

  • dx (float) – Horizontal cell size in meters.

  • dy (float) – Vertical cell size in meters.

  • dz (float) – Longitudinal cell size in meters.

  • context (XfContext) – identifies the context on which the computation is executed.

Returns

Poisson solver object.

Return type

(FFTSolver3D)

solve(rho)

Solves Poisson’s equation in free space for a given charge density.

Parameters

rho (float64 array) – charge density at the grid points in Coulomb/m^3.

Returns

electric potential at the grid points in Volts.

Return type

phi (float64 array)