Modules¶
OSCARS consists of several modules. The functions of each are found below. The sr and th modules are responsible for the vast majority of the computations. The th module is largely a set of calculations from theory whereas the sr module handles arbitrary field and beam conditions. The other modules are for plotting and simulating parametric surfaces. The OSCARS modules are:
oscars.sr¶
Primer for the Examples
All examples here assume that you have the OSCARS SR and TH modules in your path and have an OSCARS SR (TH) object called osr (oth), perhaps as in the following example
# Import the OSCARS SR module
import oscars.sr
# Creat an OSCARS SR object
osr = oscars.sr.sr()
OSCARSSR module extension.
-
oscars.sr.
add_bfield_file
([ifile, bifile, iformat, rotations, translation, scale, name]) Add a magnetic field from a text file ifile according to the format iformat.
- Currently OSCARS accepts the following file formats for iformat
- ‘OSCARS’
- ‘OSCARS1D [plus format string]’
- ‘SPECTRA’
- ‘SRW’
For the OSCARS1D format you must also include the order of the data columns, which would typically look something like: ‘OSCARS1D Z Bx By Bz’. You may use X or Y instead of Z and the order and number of B[xyz] does not matter. This mode also accepts non-uniformly distributed data.
Optionally you can rotate and translate the field in space. You can use an input scale list to scale the input (which must be in SI units of [m] for distances/positions and [T] for magnetic field values.
The rotation is performed first in the order: \(\theta_x, \theta_y, \theta_z\)
scale is a list of less than or equal length to the number of elements in iformat when OSCARS1D is selected. This will scale the input values of the i-th column before any rotation or translation. This is useful if your data file is not in [T] and [m]
Parameters: ifile : str
Name of input file
bifile : str
Name of binary input file
iformat : str
Format of the input file. Must be specified with ifile, but not with bifile
rotations : list, optional
3-element list representing rotations around x, y, and z axes: [\(\theta_x, \theta_y, \theta_z\)]
translation : list, optional
3-element list representing a translation in space [x, y, z]
scale : list, optional
List of scale factors to be used for multiplying the inputs in the order of iformat (equal in length or less than the number of input parameters)
name : str
Name of this magnetic field
Returns: None
Examples
Add a magnetic field from a data file where the data is in columns in the order Z Bx By Bz
>>> osr.add_bfield_file(ifile='file.txt', iformat='OSCARS1D Z Bx By Bz')
-
oscars.sr.
add_bfield_function
(func[, name]) Adds field in the form of a user defined python function. The input for this function must be (x, y, z, t) and return [Fx, Fy, Fz]
Parameters: function : func
A python function with input [x, y, z, t] and return of [Fx, Fy, Fz]
name : str
Name of this field
Returns: None
Examples
Create a function in python and use it as a field in OSCARS
>>> def myfunc(x, y, z, t): ... "Do not forget to write a docstring" ... if z > 0: ... return 1 ... return 0 >>> osr.add_bfield_function(myfunc)
-
oscars.sr.
add_bfield_gaussian
(bfield, sigma[, rotations, translation, name]) Add a gaussian field in 3D with the peak field magnitude and direction given by bfield, centered about a point with a given sigma in each coordinate. If any component of sigma is less than or equal to zero it is ignored (ie. spatial extent is infinite).
The functional form for this is \(\exp(-(x - x_0)^2 / \sigma_x^2)\)
Parameters: bfield : list
A list representing the peak field [Fx, Fy, Fz]
sigma : list
A list representing the sigma in each coordinate \([\sigma_x, \sigma_y, \sigma_z]\)
rotations : list, optional
3-element list representing rotations around x, y, and z axes: [\(\theta_x, \theta_y, \theta_z\)]
translation : list, optional
3-element list representing a translation in space [x, y, z]
name : str
Name of this field
Returns: None
Examples
Add a field of 0.5 [T] in the X-direction centered at X=0, Y=0, Z=0 [m], with a sigma of 0.1 [m] in the Z-direction
>>> osr.add_bfield_gaussian(bfield=[1, 0, 0], sigma=[0, 0, 0.10])
Add a field of 1 [T] in the Y-direction centered at X=1, Y=1, Z=1 [m], with a sigma of 0.05 in the Z-direction
>>> osr.add_bfield_gaussian(bfield=[0, 1, 0], sigma=[0, 0, 0.05], translation=[1, 1, 1])
-
oscars.sr.
add_bfield_interpolated
(mapping, iformat, parameter[, rotations, translation, scale, name]) Add field given a paramater and a mapping between known parameters and field data files where the field is interpolated from the known data points.
Parameters: mapping : list [[float, str], [float, str], ...]
List of parameters and associated filenames [[p1, file1], [p2, file2], ...]
iformat : str
Which input format to use (see add_bfield_file() for formats)
parameter : float
Value of parameter you are interested in
rotations : list, optional
3-element list representing rotations around x, y, and z axes: [\(\theta_x, \theta_y, \theta_z\)]
translation : list, optional
3-element list representing a translation in space [x, y, z]
scale : list, optional
List of scale factors to be used for multiplying the inputs in the order of iformat (equal in length or less than the number of input parameters)
name : str
Name of this magnetic field
Returns: None
Examples
Undulator gap interpolation between known magnetic measurement points
>>> file_list = [ ... [10.9, 'file_10.9.dat'], ... [11.0, 'file_11.0.dat'], ... [13.2, 'file_13.2.dat'], ... [16.9, 'file_16.9.dat'] ... ] >>> osr.add_bfield_interpolated(mapping=file_list, iformat='OSCARS1D Z Bx By Bz', parameter=12.123)
-
oscars.sr.
add_bfield_quadrupole
(K, width[, rotations, translation, name]) Adds a quadrupole field in a given volume according to:
\[ \begin{align}\begin{aligned}Bx = K * y\\By = K * x\\Bz = 0\end{aligned}\end{align} \]For other directions or rotated skew quad fields one can use the rotations parameter
Parameters: K : float
Quadrupole strength given in [T/m]
width : float
length in z direction (quadrupole axis direction)
rotations : list, optional
3-element list representing rotations around x, y, and z axes: [\(\theta_x, \theta_y, \theta_z\)]
translation : list, optional
3-element list representing a translation in space [x, y, z]
name : str
Name of this field
Returns: None
-
oscars.sr.
add_bfield_undulator
(bfield, period, nperiods[, phase, rotations, translation, taper, name]) Adds an ideal sinusoidal undulator field with a given maximum bfield amplitude, period, and number of periods. Optionally one can specify the phase offset (in [rad]), rotations and translation. The number of periods given is the full number of fields not counting the terminating fields.
Parameters: bfield : list
A list representing the peak field [Bx, By, Bz] in [T]
period : list
Length of one period
nperiods : int
Number of periods
phase : float
Phase offset in [rad]
rotations : list, optional
3-element list representing rotations around x, y, and z axes: [\(\theta_x, \theta_y, \theta_z\)]
translation : list, optional
3-element list representing a translation in space [x, y, z]
taper : float
The fractional change of the bfield per meter along the magnetic axis
name : str
Name of this field
Returns: None
Examples
Add an idealized undulator with 41 periods having a period of 0.050 [m] with a maximum field of 1 [T] in the y-direction where the magnetic axis is along the z-axis
>>> osr.add_bfield_undulator(bfield=[0, 1, 0], period=[0, 0, 0.050], nperiods=41)
-
oscars.sr.
add_bfield_uniform
(bfield[, width, rotations, translation, name]) Add a uniform field in a given range or for all space. The bfield is given as a 3D vector representing the field magnitude and direction. width is an optional parameters, if not present the field permeates all space. If a component of the 3D list width is less than or equal to zero, that coordinate will be ignored when calculating the field.
Parameters: bfield : list
A list representing the field [Fx, Fy, Fz]
width : list
A list representing the spetial extent of the field [Wx, Wy, Wz]
rotations : list, optional
3-element list representing rotations around x, y, and z axes: [\(\theta_x, \theta_y, \theta_z\)]
translation : list, optional
3-element list representing a translation in space [x, y, z]
name : str
Name of this field
Returns: None
Examples
Add a field of 0.0001 [T] in the X-direction for all space
>>> osr.add_bfield_uniform(bfield=[0.0001, 0, 0])
Add a field of 0.0005 [T] in the Y-direction with a width in the Z-direction of 1.5 [m] (the other directions are ignored) centered at X=0, Y=0, Z=0.75 [m].
>>> osr.add_bfield_uniform(bfield=[0, 1, 0], width=[0, 0, 1.5], translation=[0, 0, 0.75])
Add a field of 1 [T] in the X-direction in a volume of dimensions 1m x 1m x 1m.
>>> osr.add_bfield_uniform(bfield=[1, 0, 0], width=[1, 1, 1])
-
oscars.sr.
add_drift_box
(width[, rotations, translation, name]) Add a drift volume box. The particle will drift with constant velocity inside of this volume and the points inside of the volume will not be recorded in the trajectory. If a dimension of the box is <= 0 that axis is ignored.
Parameters: width: list
Width of the drift box in x, y, and z
rotations : list, optional
3-element list representing rotations around x, y, and z axes: [\(\theta_x, \theta_y, \theta_z\)]
translation : list, optional
3-element list representing a translation in space [x, y, z]
name : str
Name of the drift volume
Returns: None
-
oscars.sr.
add_efield_file
(ifile, iformat[, rotations, translation, scale, name]) Add a electric field from a text file ifile according to the format iformat.
- Currently OSCARS accepts the following file formats for iformat
- ‘OSCARS’
- ‘OSCARS1D [plus format string]’
- ‘SPECTRA’
- ‘SRW’
For the OSCARS1D format you must also include the order of the data columns, which would typically look something like: ‘OSCARS1D Z Ex Ey Ez’. You may use X or Y instead of Z and the order and number of E[xyz] does not matter. This mode also accepts non-uniformly distributed data.
Optionally you can rotate and translate the field in space. You can use an input scale list to scale the input (which must be in SI units of [m] for distances/positions and [V/m] for electric field values.
The rotation is performed first in the order: \(\theta_x, \theta_y, \theta_z\)
scale is a list of less than or equal length to the number of elements in iformat when OSCARS1D is selected. This will scale the input values of the i-th column before any rotation or translation. This is useful if your data file is not in [V/m] and [m]
Parameters: ifile : str
Name of input file
iformat : str
Format of the input file
rotations : list, optional
3-element list representing rotations around x, y, and z axes: [\(\theta_x, \theta_y, \theta_z\)]
translation : list, optional
3-element list representing a translation in space [x, y, z]
scale : list, optional
List of scale factors to be used for multiplying the inputs in the order of iformat (equal in length or less than the number of input parameters)
name : str
Name of this field
Returns: None
Examples
Add a field from a data file where the data is in columns in the order Z Ex Ey Ez
>>> osr.add_bfield_file(ifile='file.txt', iformat='OSCARS1D Z Bx By Bz')
-
oscars.sr.
add_efield_function
(func[, name]) Adds field in the form of a user defined python function. The input for this function must be (x, y, z, t) and return [Fx, Fy, Fz]
Parameters: function : func
A python function with input [x, y, z, t] and return of [Fx, Fy, Fz]
name : str
Name of this field
Returns: None
Examples
Create a function in python and use it as a field in OSCARS
>>> def myfunc(x, y, z, t): ... "Do not forget to write a docstring" ... if z > 0: ... return 1 ... return 0 >>> osr.add_efield_function(myfunc)
-
oscars.sr.
add_efield_gaussian
(efield, sigma[, rotations, translation]) Add a gaussian field in 3D with the peak field magnitude and direction given by efield, centered about a point with a given sigma in each coordinate. If any component of sigma is less than or equal to zero it is ignored (ie. spatial extent is infinite).
The functional form for this is \(\exp(-(x - x_0)^2 / \sigma_x^2)\)
Parameters: efield : list
A list representing the peak field [Fx, Fy, Fz]
sigma : list
A list representing the sigma in each coordinate \([\sigma_x, \sigma_y, \sigma_z]\)
rotations : list, optional
3-element list representing rotations around x, y, and z axes: [\(\theta_x, \theta_y, \theta_z\)]
translation : list, optional
3-element list representing a translation in space [x, y, z]
Returns: None
Examples
Add a field of 0.5 [V/m] in the X-direction centered at X=0, Y=0, Z=0 [m], with a sigma of 0.1 [m] in the Z-direction
>>> osr.add_efield_gaussian(efield=[1, 0, 0], sigma=[0, 0, 0.10])
Add a field of 1 [V/m] in the Y-direction centered at X=1, Y=1, Z=1 [m], with a sigma of 0.05 in the Z-direction
>>> osr.add_efield_gaussian(efield=[0, 1, 0], sigma=[0, 0, 0.05], translation=[1, 1, 1])
-
oscars.sr.
add_efield_interpolated
(mapping, iformat, parameter[, rotations, translation, scale, name]) Add field given a paramater and a mapping between known parameters and field data files where the field is interpolated from the known data points.
Parameters: mapping : list [[float, str], [float, str], ...]
List of parameters and associated filenames [[p1, file1], [p2, file2], ...]
iformat : str
Which input format to use (see add_efield_file() for formats)
parameter : float
Value of parameter you are interested in
rotations : list, optional
3-element list representing rotations around x, y, and z axes: [\(\theta_x, \theta_y, \theta_z\)]
translation : list, optional
3-element list representing a translation in space [x, y, z]
scale : list, optional
List of scale factors to be used for multiplying the inputs in the order of iformat (equal in length or less than the number of input parameters)
name : str
Name of this field
Returns: None
Examples
See add_bfield_interpolated() for examples
-
oscars.sr.
add_efield_undulator
(efield, period, nperiods[, phase, rotations, translation, taper]) Adds an ideal sinusoidal undulator field with a given maximum efield amplitude, period, and number of periods. Optionally one can specify the phase offset (in [rad]), rotations and translation. The number of periods given is the full number of fields not counting the terminating fields.
Parameters: efield : list
A list representing the peak field [Fx, Fy, Fz] in [V/m]
period : list
Length of one period
nperiods : int
Number of periods
phase : float
Phase offset in [rad]
rotations : list, optional
3-element list representing rotations around x, y, and z axes: [\(\theta_x, \theta_y, \theta_z\)]
translation : list, optional
3-element list representing a translation in space [x, y, z]
taper : float
The fractional change of the efield per meter along the magnetic axis
name : str
Name of this field
Returns: None
Examples
Add an idealized undulator with 41 periods having a period of 0.050 [m] with a maximum field of 1 [T] in the y-direction where the magnetic axis is along the z-axis
>>> osr.add_efield_undulator(efield=[0, 1, 0], period=[0, 0, 0.050], nperiods=41)
-
oscars.sr.
add_efield_uniform
(efield[, width, rotations, translation, name]) Add a uniform field in a given range or for all space. The efield is given as a 3D vector representing the field magnitude and direction. width is an optional parameters, if not present the field permeates all space. If a component of the 3D list width is less than or equal to zero, that coordinate will be ignored when calculating the field.
Parameters: efield : list
A list representing the field [Fx, Fy, Fz]
width : list
A list representing the spetial extent of the field [Wx, Wy, Wz]
rotations : list, optional
3-element list representing rotations around x, y, and z axes: [\(\theta_x, \theta_y, \theta_z\)]
translation : list, optional
3-element list representing a translation in space [x, y, z]
name : str
Name of this field
Returns: None
Examples
Add a field of 0.0001 [V/m] in the X-direction for all space
>>> osr.add_efield_uniform(efield=[0.0001, 0, 0])
Add a field of 0.0005 [V/m] in the Y-direction with a width in the Z-direction of 1.5 [m] (the other directions are ignored) centered at X=0, Y=0, Z=0.75 [m].
>>> osr.add_efield_uniform(efield=[0, 1, 0], width=[0, 0, 1.5], translation=[0, 0, 0.75])
Add a field of 1 [V/m] in the X-direction in a volume of dimensions 1m x 1m x 1m.
>>> osr.add_efield_uniform(efield=[1, 0, 0], width=[1, 1, 1])
-
oscars.sr.
add_particle_beam
([type, name, energy_GeV, d0, x0, beam, sigma_energy_GeV, t0, current, weight, rotations, translation, horizontal_direction, beta, alpha, gamma, emittance, eta, lattice_reference, mass, charge]) Add a particle beam to the OSCARS object with a name given by name. There is no limit to the number of different particle beams one can add. They are added with a weight which is by default 1. The weight is used in random sampling when asking for a new particle, for example in oscars.sr.set_new_particle(). If the beam parameter is given you only need to specify name and x0.
- Supported particle types for type are:
- electron
- positron
- muon
- anti-muon
- proton
- anti-proton
- pi+
- pi-
Parameters: type : str
One of the built-in types of particle beams, or ‘custom’. If you use custom you must also specify mass and charge.
name : str
User identified of this beam
energy_GeV : float
Beam energy in [GeV]
d0 : list
Vector [float, float, float] representing the default direction of the beam at the initial position. The normalization of this vector does not matter.
x0 : list
Coordinates of the initial position in [m] [x, y, z]
beam : str
Name of predefined beam
sigma_energy_GeV : float
Beam energy spread in [GeV]
t0 : float
Initial time in [m] at the initial_position. Time here is in terms of ct.
current : float
Beam current in [A]. If this parameter is 0, the current defaults to the single particle charge
weight : float
Weight to give this beam for random sampling when there are multiple beams
rotations : list, optional
3-element list representing rotations around x, y, and z axes: [\(\theta_x, \theta_y, \theta_z\)]
translation : list, optional
3-element list representing a translation in space [x, y, z]
horizontal_direction : list
A list representing the horizontal beam direction. This must be perpendicular to the beam direction. The vertical direction is defined internally as \([direction] \times [horizontal\_direction]\).
beta : list
Values of the horizontal and vertical beta funtion at the point lattice_center [beta_x, beta_y]
alpha : list
Values of the horizontal and vertical alpha funtion at the point lattice_center [alpha_x, alpha_y]
gamma : list
Values of the horizontal and vertical gamma funtion at the point lattice_center [gamma_x, gamma_y]
emittance : list
values of the horizontal and vertical emittance [emittance_x, emittance_y]
eta : list
Values of the horizontal and vertical dispersion at the point lattice_center [eta_x, eta_y]. Currently eta_y is ignored and dispersion is only used in oscars.th calculations
lattice_reference : list
Coordinates of the lattice center [x, y, z] (must be on-axis with respect to the beam)
mass : float
mass of a custom particle. This is only used if type = ‘custom’. Must be non-zero.
charge : float
Charge of a custom particle. This is only used if type = ‘custom’. Must be non-zero.
Returns: None
Currently the predefined beam types are:
- NSLSII
- NSLSII-ShortStraight
- NSLSII-LongStraight
Examples
Add an electron beam with 0.500 [A] current at an initial position of [0, 0, 0] in the Y direction with energy of 3 [GeV]
>>> osr.add_particle_beam(type='electron', name='beam_0', x0=[0, 0, 0], d0=[0, 1, 0], energy_GeV=3, current=0.500)
Add a positron beam with 0.500 [A] current at an initial position of [-2, 0, 0] in the direction given by theta in the X-Y plane with energy of 3 [GeV]
>>> from math import sin, cos >>> theta = 0.25 * osr.pi() >>> osr.add_particle_beam(type='positron', name='beam_0', x0=[-2, 0, 0], d0=[sin(theta), cos(theta), 0], energy_GeV=3, current=0.500)
Add a predefined beam for the NSLSII short straight section
>>> osr.add_particle_beam(beam='NSLSII-ShortStraight', name='beam_0', x0=[-2, 0, 0])
-
oscars.sr.
add_to_flux
(flux[, weight]) Add flux map to the current flux map with a given weight.
Parameters: flux : list
A list of points and fluxes (in flux format)
weight : float
Weight for this flux map
Returns: None
-
oscars.sr.
add_to_power_density
(power_density[, weight]) Add power_density map to the current power_density map with a given weight.
Parameters: power_density : list
A list of points and power_densities (in power_density format)
weight : float
Weight for this power_density map
Returns: None
-
oscars.sr.
add_to_spectrum
(spectrum[, weight]) Add a spectrum to the current spectrum with a given weight.
Parameters: spectrum : list
A list of pairs of numbers (spectrum format)
weight : float
Weight for this spectrum
Returns: None
-
oscars.sr.
average_flux
([ifiles, bifiles, ofile, bofile, dim]) Average from different files and output to specified file. The input files must have the same format.
Parameters: ifiles : list
The input file names as strings: [‘f0.txt’, ‘f1.txt’, ...]
bifiles : list
The binary input file names as strings: [‘f0.dat’, ‘f1.dat’, ...]
ofile : str
The output file name
bofile : str
The binary output file name
dim : int
in 2 or 3 dimensions (default is 2)
Returns: flux : list
A list, each element of which is a pair representing the position (x) and value (v) at that position. eg [[[x1_0, x2_0, x3_0], v_0], [[x1_1, x2_1, x3_1], v_1]], ...]. The position is always given as a list of length 3.
-
oscars.sr.
average_power_density
([ifiles, bifiles, ofile, bofile, dim]) Average from different files and output to specified file. The input files must have the same format.
Parameters: ifiles : list
The input file names as strings: [‘f0.txt’, ‘f1.txt’, ...]
bifiles : list
The binary input file names as strings: [‘f0.dat’, ‘f1.dat’, ...]
ofile : str
The output file name
bofile : str
The binary output file name
dim : int
in 2 or 3 dimensions (default is 2)
Returns: power_density : list
A list, each element of which is a pair representing the position (x) and value (v) at that position. eg [[[x1_0, x2_0, x3_0], v_0], [[x1_1, x2_1, x3_1], v_1]], ...]. The position is always given as a list of length 3.
-
oscars.sr.
average_spectra
([ifiles, bifiles, ofile, bofile]) Average spectra from different files. The input files must have the same format.
Parameters: ifiles : list
The input file names as strings: [‘f0.txt’, ‘f1.txt’, ...]
bifiles : list
The binary input file names as strings: [‘f0.dat’, ‘f1.dat’, ...]
ofile : str
The output file name
bofile : str
The binary output file name
Returns: spectrum : list
A list of 2D lists, each of which is a pair representing the energy and flux at that energy. eg [[energy_0, flux_0], [energy_1, flux_1], ...]
-
oscars.sr.
calculate_efield_vs_time
(obs[, ofile]) Calculate the electric field in the time domain for a single particle
See the Mathematical Notes section for the expression used in this calculation.
Parameters: obs : lsit
Point where you wish to calculate the electric field [x, y, z]
ofile : str
Output file name
Returns: efield : list
A list, each element of which has a time (in [s]) and a 3-dimensional list representing the x, y, and z componemts of the electric field at that time: [[t, [Ex, Ey, Ez]], ...]
-
oscars.sr.
calculate_flux
(energy_eV, points[, normal, rotations, translation, nparticles, nthreads, gpu, ngpu, precision, max_level, max_level_extended, ofile, bofile, quantity]) Calculates the flux at a given set of points
See the Mathematical Notes section for the expression used in this calculation.
Parameters: energy_eV : float
Photon energy of interest
points : list
A list of points, each point containing a position in 3D (as a list) and a normal vector at that position (also as a 3D list): [[[x, y, z], [nx. ny. nz]], [...], ...]
normal : int
-1 if you wish to reverse the normal vector, 0 if you wish to ignore the +/- direction in computations, 1 if you with to use the direction of the normal vector as given.
rotations : list, optional
3-element list representing rotations around x, y, and z axes: [\(\theta_x, \theta_y, \theta_z\)]
translation : list, optional
3-element list representing a translation in space [x, y, z]
nparticles : int
Number of particles to use for multi-particle calculations
nthreads : int
Number of threads to use
gpu : int
Use the gpu or not (0 or 1). If 1 will attempt to use ALL gpus available. This is overridden if you use the input ‘ngpu’
ngpu : int or list
If ngpu is an int, use that number of gpus (if available). If ngpu is a list, the list should be a list of gpus you wish to use
precision : float
Calculation precision parameter (typically 0.01 which is 1%)
max_level: int
Maximum “level” to us for trajectory in the calculation. Level N corresponds to a total of 2**(N+2) trajectory points. You cannot go beyond the internal maximum. You are not guaranteed precision parameter is met if this is used.
max_level_extended: int
Maximum “level” to us for trajectory in the calculation. If set to higher than max_level the computation will proceed beyond max_level without creating trajectory arrays in memory (but it will be slower)
ofile : str
Output file name
bofile : str
Binary output file name
quantity: str
Quantity to return. Available are:
‘power density’ (default) ‘precision’ - Estimated precision for each point ‘level’ - Trajectory level reached (npoints = 2**(n+1) - 1), if return is -1 the requested precision was not reached
Returns: power_density : list
A list, each element of which is a pair representing the position (2D relative (default) or 3D absolute) and power density [\(W / mm^2\)] at that position. eg [[[x1_0, x2_0, x3_0], pd_0], [[x1_1, x2_1, x3_1], pd_1]], ...]. The position is always given as a list of length 3. For the default (dim=2) the third element is always zero.
-
oscars.sr.
calculate_flux_rectangle
(energy_eV, npoints[, plane, normal, dim, width, rotations, translation, x0x1x2, polarization, angle, horizontal_direction, propogation_direction, nparticles, nthreads, gpu, ngpu, precision, max_level, max_level_extended, quantity, ofile, bofile]) Calculate the flux density in a rectangle either defined by three points, or by defining the plane the rectangle is in and the width, and then rotating and translating it to where it needs be. The simplest is outlined in the first example below. By default (dim=2) this returns a list whose position coordinates are in the local coordinate space x1 and x2 (ie they do not include the rotations and translation). if dim=3 the coordinates in the return list are in absolute 3D space.
You must specify either both (plane and width) or x0x1x2
See the Mathematical Notes section for the expression used in this calculation.
Parameters: energy_eV : float
Photon energy of interest
npoints : list [int, int]
Number of points in X1 and X2 dimension [n1, n2]
plane : str
The plane to start in (XY, XZ, YZ, YX, ZX, ZY). The normal to the surface is defined using the right handed cross product (ie the last three have opposite normal vectors from the first three)
normal : int
-1 if you wish to reverse the normal vector, 0 if you wish to ignore the +/- direction in computations, 1 if you with to use the direction of the normal vector as given.
dim : int
Defaults to 2 where output is in the local plane coordinates X1 and X2. If you want the return to be given in 3D set dim=3 which will return with X, Y, and Z in absolute coordinates.
width : list
Width of rectangle in X1 and X2: [w1, w2]
rotations : list, optional
3-element list representing rotations around x, y, and z axes: [\(\theta_x, \theta_y, \theta_z\)]
translation : list, optional
3-element list representing a translation in space [x, y, z]
x0x1x2 : list
List of three points [[x0, y0, z0], [x1, y1, z1], [x2, y2, z2]] defining a parallelogram (vectors 0->1, and 0->2)
polarization : str
Which polarization mode to calculate. Can be ‘all’, ‘linear-horizontal’, ‘linear-vertical’, ‘circular-left’, ‘circular-right’, or ‘linear’ (if linear you must specify the angle parameter)
angle : float
Only used if polarization=’linear’ is specified. The ‘angle’ is that from the horizontal_direction for the polarization directino you are interested in
horizontal_direction : list
The direction you consider to be horizontal. Should be perpendicular to the photon beam propogation direction
propogation_direction : list
Propogation direction of photon beam
nparticles : int
The number of particles you wish to run for a multi-particle simulation
nthreads : int
Number of threads to use
gpu : int
Use the gpu or not (0 or 1). If 1 will attempt to use ALL gpus available. This is overridden if you use the input ‘ngpu’
ngpu : int or list
If ngpu is an int, use that number of gpus (if available). If ngpu is a list, the list should be a list of gpus you wish to use
precision : float
Calculation precision parameter (typically 0.01 which is 1%)
max_level: int
Maximum “level” to use for trajectory in the calculation. Level N corresponds to a total of 2**(N+2) trajectory points. You cannot go beyond the internal maximum. You are not guaranteed precision parameter is met if this is used.
max_level_extended: int
Maximum “level” to use for trajectory in the calculation. If set to higher than max_level the computation will proceed beyond max_level without creating trajectory arrays in memory (but it will be slower)
quantity: str
Quantity to return. Available are:
‘flux’ (default) ‘precision’ - Estimated precision for each point ‘level’ - Trajectory level reached (npoints = 2**(n+1) - 1), if return is -1 the requested precision was not reached
ofile : str
Output file name
bofile : str
Binary output file name
Returns: flux : list
A list, each element of which is a pair representing the position (2D relative (default) or 3D absolute) and flux [\(W / mm^2\)] at that position. eg [[[x1_0, x2_0, x3_0], f_0], [[x1_1, x2_1, x3_1], f_1]], ...]. The position is always given as a list of length 3. For the default (dim=2) the third element is always zero.
-
oscars.sr.
calculate_power_density
(points[, normal, rotations, translation, nparticles, gpu, nthreads, precision, max_level, max_level_extended, quantity, ofile]) Calculate the power density for each point in the list points.
See the Mathematical Notes section for the expression used in this calculation.
Parameters: points : list
A list of points, each point containing a position in 3D (as a list) and a normal vector at that position (also as a 3D list): [[[x, y, z], [nx. ny. nz]], [...], ...]
normal : int
-1 if you wish to reverse the normal vector, 0 if you wish to ignore the +/- direction in computations, 1 if you with to use the direction of the normal vector as given.
rotations : list, optional
3-element list representing rotations around x, y, and z axes: [\(\theta_x, \theta_y, \theta_z\)]
translation : list, optional
3-element list representing a translation in space [x, y, z]
nparticles : int
Number of particles to use for multi-particle calculations
gpu : int
Use the gpu or not (0 or 1)
nthreads : int
Number of threads to use
precision : float
Calculation precision parameter (typically 0.01 which is 1%)
max_level: int
Maximum “level” to use for trajectory in the calculation. Level N corresponds to a total of 2**(N+2) trajectory points. You cannot go beyond the internal maximum. You are not guaranteed precision parameter is met if this is used.
max_level_extended: int
Maximum “level” to use for trajectory in the calculation. If set to higher than max_level the computation will proceed beyond max_level without creating trajectory arrays in memory (but it will be slower)
quantity: str
Quantity to return. Available are:
‘power density’ (default) ‘precision’ - Estimated precision for each point ‘level’ - Trajectory level reached (npoints = 2**(n+1) - 1), if return is -1 the requested precision was not reached
ofile : str
Output file name
Returns: power_density : list
A list, each element of which is a pair representing the position (2D relative (default) or 3D absolute) and power density [\(W / mm^2\)] at that position. eg [[[x1_0, x2_0, x3_0], pd_0], [[x1_1, x2_1, x3_1], pd_1]], ...]. The position is always given as a list of length 3. For the default (dim=2) the third element is always zero.
-
oscars.sr.
calculate_power_density_line
(x1, x2, normal_direction, npoints[, ofile, bofile, normal, nparticles, gpu, nthreads, precision, max_level]) See the Mathematical Notes section for the expression used in this calculation.
Parameters: x1 : list[3]
Point in space to start [x0, y0, z0]
x2 : list[3]
Point in space to end [x1, y1, z1]
normal_direction : list[3]
Direction of the normal vector for power density calculation. Typically this is
npoints: int
number of in each dimension for surface
ofile : str
Output file name
bofile : str
Binary output file name
normal : int
-1 if you wish to reverse the normal vector, 0 if you wish to ignore the +/- direction in computations, 1 if you with to use the direction of the normal vector as given.
nparticles : int
Number of particles to use for multi-particle calculations
gpu : int
Use the gpu or not (0 or 1)
nthreads : int
Number of threads to use
precision : float
Calculation precision parameter (typically 0.01 which is 1%)
max_level: int
Maximum “level” to us for trajectory in the calculation. Level N corresponds to a total of 2**(N+2) trajectory points. You cannot go beyond the internal maximum. You are not guaranteed precision parameter is met if this is used.
max_level_extended: int
Maximum “level” to use for trajectory in the calculation. If set to higher than max_level the computation will proceed beyond max_level without creating trajectory arrays in memory (but it will be slower)
Returns: power_density_1d : list
A list, each element of which is a pair representing the length along the line specified and power density [\(W / mm^2\)] at that position. eg [[x0, p0], [x1, p1], ...]
-
oscars.sr.
calculate_power_density_rectangle
(npoints[, plane, width, x0x1x2, rotations, translation, ofile, bofile, normal, nparticles, gpu, ngpu, nthreads, precision, max_level, max_level_extended, dim, quantity]) Calculate the power density in a rectangle either defined by three points, or by defining the plane the rectangle is in and the width, and then rotating and translating it to where it needs be. The simplest is outlined in the first example below. By default (dim=2) this returns a list whose position coordinates are in the local coordinate space x1 and x2 (ie they do not include the rotations and translation). if dim=3 the coordinates in the return list are in absolute 3D space.
See the Mathematical Notes section for the expression used in this calculation.
You must specify either both (plane and width) or x0x1x2
Parameters: npoints: int
number of in each dimension for surface
plane : str
The plane to start in (XY, XZ, YZ, YX, ZX, ZY). The normal to the surface is defined using the right handed cross product (ie the last three have opposite normal vectors from the first three)
width : list
Width of rectangle in X1 and X2: [w1, w2]
x0x1x2 : list
List of three points [[x0, y0, z0], [x1, y1, z1], [x2, y2, z2]] defining a parallelogram (vectors 0->1, and 0->2)
rotations : list, optional
3-element list representing rotations around x, y, and z axes: [\(\theta_x, \theta_y, \theta_z\)]
translation : list, optional
3-element list representing a translation in space [x, y, z]
ofile : str
Output file name
bofile : str
Binary output file name
normal : int
-1 if you wish to reverse the normal vector, 0 if you wish to ignore the +/- direction in computations, 1 if you with to use the direction of the normal vector as given.
nparticles : int
Number of particles to use for multi-particle calculations
gpu : int
Use the gpu or not (0 or 1). If 1 will attempt to use ALL gpus available. This is overridden if you use the input ‘ngpu’
ngpu : int or list
If ngpu is an int, use that number of gpus (if available). If ngpu is a list, the list should be a list of gpus you wish to use
nthreads : int
Number of threads to use
precision : float
Calculation precision parameter (typically 0.01 which is 1%)
max_level: int
Maximum “level” to use for trajectory in the calculation. Level N corresponds to a total of 2**(N+2) trajectory points. You cannot go beyond the internal maximum. You are not guaranteed precision parameter is met if this is used.
max_level_extended: int
Maximum “level” to use for trajectory in the calculation. If set to higher than max_level the computation will proceed beyond max_level without creating trajectory arrays in memory (but it will be slower)
dim : int
Defaults to 2 where output is in the local plane coordinates X1 and X2. If you want the return to be given in 3D set dim=3 which will return with X, Y, and Z in absolute coordinates.
quantity: str
Quantity to return. Available are:
‘power density’ (default) ‘precision’ - Estimated precision for each point ‘level’ - Trajectory level reached (npoints = 2**(n+1) - 1), if return is -1 the requested precision was not reached
Returns: power_density : list
A list, each element of which is a pair representing the position (2D relative or 3D absolute) and power density [\(W / mm^2\)] at that position. eg [[x1_0, x2_0], pd_0, [x1_1, x2_1], pd_1], ...]
Examples
Calculate the power density within a simple rectangle 1 [cm] x 1 [cm], 30 [m] downstream
>>> osr.calculate_power_density(plane='XY', width=[0.01, 0.01], npoints=[51, 51], translation=[0, 0, 30])
Calculate the power density within a simple rectangle 1 [cm] x 1 [cm], 30 [m] downstream defined using the three-point method
>>> osr.calculate_power_density(x0x1x2=[[-0.005, -0.005, 30], [+0.005, -0.005, 30], [-0.005, +0.005, 30]], npoints=[51, 51])
Calculate the power density on a flat surface close to and parallel to the beam direction (beam in z-direction)
>>> power_density = osr.calculate_power_density(plane='XZ', width=[0.008, 2], npoints=[51, 101], normal=-1)
-
oscars.sr.
calculate_power_density_stl
(ifiles[, rotations, translation, ofile, bofile, stlofile, normal, scale, nparticles, gpu, ngpu, nthreads, precision, max_level, max_level_extended, quantity]) Calculate the power density on surfaces described in STL format input files
See the Mathematical Notes section for the expression used in this calculation.
You must specify either both (plane and width) or x0x1x2
Parameters: ifiles : list[str]
A list of STL files to import.
rotations : list, optional
3-element list representing rotations around x, y, and z axes: [\(\theta_x, \theta_y, \theta_z\)]
translation : list, optional
3-element list representing a translation in space [x, y, z]
ofile : str
Output file name
bofile : str
Binary output file name
stlofile : str
STL output file name. Will include RGB color information in the attribute byte count, maybe
normal : int
-1 if you wish to reverse the normal vector, 0 if you wish to ignore the +/- direction in computations, 1 if you with to use the direction of the normal vector as given.
scale : float
What to scale the dimensions by. Default input is in meters.
nparticles : int
Number of particles to use for multi-particle calculations
gpu : int
Use the gpu or not (0 or 1). If 1 will attempt to use ALL gpus available. This is overridden if you use the input ‘ngpu’
ngpu : int or list
If ngpu is an int, use that number of gpus (if available). If ngpu is a list, the list should be a list of gpus you wish to use
nthreads : int
Number of threads to use
precision : float
Calculation precision parameter (typically 0.01 which is 1%)
max_level: int
Maximum “level” to use for trajectory in the calculation. Level N corresponds to a total of 2**(N+2) trajectory points. You cannot go beyond the internal maximum. You are not guaranteed precision parameter is met if this is used.
max_level_extended: int
Maximum “level” to use for trajectory in the calculation. If set to higher than max_level the computation will proceed beyond max_level without creating trajectory arrays in memory (but it will be slower)
quantity: str
Quantity to return. Available are:
‘power density’ (default) ‘precision’ - Estimated precision for each point ‘level’ - Trajectory level reached (npoints = 2**(n+1) - 1), if return is -1 the requested precision was not reached
Returns: power_density : list
This return list is NOT the same as other power density calculations. The return format is as follows: [ [[[T0x, T0y, T0z], [T1x, T1y, T1z], [T2x, T2y, T2z]], PD], [...], ...] where the Ts are the triangle vertices, and PD is the power density in the center of that triangle.
Examples
coming
-
oscars.sr.
calculate_spectrum
(obs[, npoints, energy_range_eV, energy_points_eV, points_eV, polarization, angle, horizontal_direction, propogation_direction, precision, max_level, nparticles, nthreads, gpu, ngpu, quantity, ofile, bofile]) Calculate the spectrum given a point in space, the range in energy, and the number of points. The calculation uses the current particle and its initial conditions. If the trajectory has not been calculated it is calculated first. The units of this calculation are [\(photons / mm^2 / 0.1% bw / s\)]
Previously to calling this function you must define a particle beam and define start and stop times at very minimum.
You must provide either (npoints and energy_range_eV) or points_eV.
Parameters: obs : list
Point [x, y, z] where you wish to calculate the spectrum
npoints : int
Number of points to calculate in the given energy range
energy_range_eV : list
energy range [min, max] in eV as a list of length 2
points_eV : list
A list of points to calculate the flux at ie [12.3, 45.6, 78.9, 123.4]
polarization : str
Which polarization mode to calculate. Can be ‘all’, ‘linear-horizontal’, ‘linear-vertical’, ‘circular-left’, ‘circular-right’, or ‘linear’ (if linear you must specify the angle parameter)
horizontal_direction : list
The direction you consider to be horizontal. Should be perpendicular to the photon beam propogation direction
vertical_direction : list
Same as horizontal_direction but the vertical direction
precision : float
Calculation precision parameter (typically 0.01 which is 1%)
max_level: int
Maximum “level” to use for trajectory in the calculation. Level N corresponds to a total of 2**(N+2) trajectory points. You cannot go beyond the internal maximum. You are not guaranteed precision parameter is met if this is used.
max_level_extended: int
Maximum “level” to use for trajectory in the calculation. If set to higher than max_level the computation will proceed beyond max_level without creating trajectory arrays in memory (but it will be slower)
angle : float
Only used if polarization=’linear’ is specified. The ‘angle’ is that from the horizontal_direction for the polarization directino you are interested in
nparticles : int
The number of particles you wish to run for a multi-particle simulation
nthreads : int
Number of threads to use
gpu : int
Use the gpu or not (0 or 1). If 1 will attempt to use ALL gpus available. This is overridden if you use the input ‘ngpu’
ngpu : int or list
If ngpu is an int, use that number of gpus (if available). If ngpu is a list, the list should be a list of gpus you wish to use
quantity: str
Quantity to return. Available are:
‘flux’ (default) ‘precision’ - Estimated precision for each point ‘level’ - Trajectory level reached (npoints = 2**(n+1) - 1), if return is -1 the requested precision was not reached
ofile : str
Output file name
bofile : str
Binary output file name
Returns: spectrum : list
A list of 2D lists, each of which is a pair representing the energy [eV] and flux [\(photons / mm^2 / 0.1% bw / s\)] at that energy. eg [[energy_0, flux_0], [energy_1, flux_1], ...]
Examples
Calculate the spectrum at a point 30 [m] downstream (assuming beam direction is in the Z direction) in the energy range 100 to 1000 [eV] with 900 points.
>>> osr.calculate_spectrum(obs=[0, 0, 30], energy_range_eV=[100, 1000], npoints=900)
-
oscars.sr.
calculate_total_power
() Calculate the total radiated power based on the current particle and that particle’s beam current.
See the Mathematical Notes section for the expression used in this calculation.
Returns: power : float
Total power in [W]
-
oscars.sr.
calculate_trajectory
() Calculates the trajectory for the current internal particle. This calculates the trajectory in 3D from the time set by oscars.sr.set_ctstart() to the time set by oscars.sr.set_ctstop() beginning at the t0 given by the particle beam from which this particle comes from. It first does a forward propogation to the stop time, then a backward propogation to the start time.
It is not necessary to call this method before other calculations such as spectrum, power density, or flux calculation methods.
If you have a current particle loaded using SetNewParticle this method will calculate the trajectory for that particle. If no particle is defined one will be randomly selected based on the beam weights and beam parameters.
Parameters: None
Returns: trajectory : list
A list of points of the form [[[x, y, z], [Beta_x, Beta_y, Beta_z]], ...]
-
oscars.sr.
check_gpu
() Will return the number of GPUs available, or a negative number for any other error, for instance if your distribution was not compiled with GPU support.
Returns: ngpu : int
-
oscars.sr.
clear_bfields
() Remove all of the existing magnetic fields
Parameters: None Returns: None
-
oscars.sr.
clear_drifts
() Remove all drift volumes
Parameters: None Returns: None
-
oscars.sr.
clear_efields
() Remove all of the existing electric fields.
Parameters: None Returns: None
-
oscars.sr.
clear_particle_beams
() Remove all of the existing particle beams
Parameters: None Returns: None
-
oscars.sr.
correct_trajectory
() Correct the trajectory to desired position and direction and specified point based on corrector kicks at specified locations.
This function will add magnetic fields to the sr object with names that begin with an underscore.
Parameters: None
Returns: trajectory : list
A list of points of the form [[[x, y, z], [Beta_x, Beta_y, Beta_z]], ...]
-
oscars.sr.
get_bfield
(x) Get the 3D field at any point in space. This is the sum of all fields added to this OSCARS object.
Parameters: x : list
A 3D list representing a point in space [x, y, z]
Returns: bfield : [float, float, float]
list representing the field [Fx, Fy, Fz]
-
oscars.sr.
get_ctstart
() Gets the start time for particle propogation. The time is in units of [m].
Returns: time, float
-
oscars.sr.
get_ctstop
() Gets the stop time for particle propogation. The time is in units of [m].
Returns: time : float
-
oscars.sr.
get_efield
(x) Get the 3D field at any point in space. This is the sum of all fields added to this OSCARS object.
Parameters: x : list
A 3D list representing a point in space [x, y, z]
Returns: bfield : [float, float, float]
list representing the field [Fx, Fy, Fz]
-
oscars.sr.
get_flux
() Get the current flux map stored in memory
Parameters: None
Returns: flux : list
list of [[[x, y, z], flux], ...]
-
oscars.sr.
get_npoints_trajectory
() Get the number of points to be used for the trajectory calculation
Returns: npoints : int
-
oscars.sr.
get_particle_beta0
() Get the initial \(\vec \beta\) for the current particle
Returns: b0 : [float, float, float]
Initial beta of particle
-
oscars.sr.
get_particle_e0
() Get the initial energy for the current particle in [GeV]
Returns: e0 : float
Initial energy of current particle
-
oscars.sr.
get_particle_x0
() Get the initial position for the current particle in [m]
Returns: x0 : [float, float, float]
Initial particle position
-
oscars.sr.
get_power_density
() Get the current power_density map stored in memory
Parameters: None
Returns: power_density : list
list of [[[x, y, z], power_density], ...]
-
oscars.sr.
get_spectrum
() Get the current spectrum
Parameters: None
Returns: spectrum : list
A list of 2D lists, each of which is a pair representing the energy and flux at that energy. eg [[energy_0, flux_0], [energy_1, flux_1], ...]
-
oscars.sr.
get_trajectory
() Get the current trajectory. If the trajectory has not been calculated this will return an empty list. The format of the returned list consists of a list of lists giving you the position and beta (v/c) of the particle at each position. For a trajectory returnned is of the form: [[[x, y, z], [Beta_x, Beta_y, Beta_z]], ...]
Parameters: None
Returns: trajectory : list
A list of points of the form [[[x, y, z], [Beta_x, Beta_y, Beta_z]], ...]
-
oscars.sr.
me
() Get the value of the mass of the electron 9.10938356e-31 [kg]
Returns: electron mass : float
-
oscars.sr.
norm
() Random float from the normal distribution centered at 0 with a sigma of 1
Returns: float
-
oscars.sr.
pi
() The value of pi
Returns: pi : float
-
oscars.sr.
print_all
() Print all internal information related to beams and fields
Parameters: None Returns: None
-
oscars.sr.
print_bfields
() Print information about all magnetic fields to standard out
Parameters: None Returns: None
-
oscars.sr.
print_drifts
() Print information about all drift volumes to the standard out
Parameters: None Returns: None
-
oscars.sr.
print_efields
() Print information about all electric fields to the standard out
Parameters: None Returns: None
-
oscars.sr.
print_gpu
() Print information about all gpus to standard out
Parameters: None Returns: None
-
oscars.sr.
print_particle_beams
() Print information for all internal particle beams
Parameters: None Returns: None
-
oscars.sr.
qe
() Get the value of elementary charge: +1.602176462e-19 [C]
Returns: charge of electron : float
-
oscars.sr.
rand
() Uniformally distributed random float in the range [0, 1)
Returns: float
-
oscars.sr.
remove_bfield
(name) Remove all fields with the given name
Parameters: name : str
Name of field to remove
Returns: None
-
oscars.sr.
remove_drift
(name) Remove all drift volumes with the given name
Parameters: name : str
Name of the drift volume
Returns: None
-
oscars.sr.
remove_efield
(name) Remove all fields with the given name
Parameters: name : str
Name of field to remove
Returns: None
-
oscars.sr.
set_ctstartstop
(start, stop) Set the start and stop times for the trajectory calculation. Start time must be less than or equal to the T0 defined in the initial conditions for the particle beam. Similarly, stop time must be greater than that T0. This also sets a default value for the number of points used in the trajectory calculation.
Parameters: start : float
Start time
stop : float
Stop time
Returns: None
-
oscars.sr.
set_gpu_global
(gpu) If set to 1, OSCARS will try to use the GPU for all calculations
Parameters: gpu : int
This will set the gpu active
Returns: None
-
oscars.sr.
set_new_particle
([beam, particle]) If no arguments are given sets the current internal particle to a random new particle. The randomization is based on the weights given for each beam. This also sets the initial conditions for the particle used in trajectory calculations based on the beam parameters within the randomly sepected beam. You can specify which beam you want a random particle from using the beam parameter. The particle parameter can be ‘ideal’ if you want the ideal initial conditions for a particle without randomization.
Parameters: beam : str
The name of the beam from which to get a particle from
particle : str
‘ideal’ or ‘random’, for random, you may omit this.
Returns: None
-
oscars.sr.
set_npoints_per_meter_trajectory
(npoints) Sets the number of points per meter to be used in the trajectory calculation.
Parameters: npoints : int
Number of points per neter
Returns: None
-
oscars.sr.
set_npoints_trajectory
(npoints) Sets the number of points to be used in the trajectory calculation. Only call this function after set_ctstartstop().
Parameters: npoints : int
Number of points
Returns: None
example:
# First set the start and stop times osr.set_ctstartstop(0, 2) # Use 12345 points in the trajectory calculations osr.set_npoints_trajectory(12345)
-
oscars.sr.
set_nthreads_global
(nthreads) Set the number of threads you wish to use for all calculations. If the GPU is requested it will take precedence over this multi-threading.
Parameters: nthreads : int
Default number of threads to use for calculations
Returns: None
-
oscars.sr.
set_particle_beam
([type, name, energy_GeV, d0, x0, beam, sigma_energy_GeV, t0, current, weight, rotations, translation, horizontal_direction, beta, alpha, gamma, emittance, eta, lattice_reference, mass, charge]) This function is the same as add_particle_beam(), but it clears all particle beams before the ‘add’.
-
oscars.sr.
set_particle_beam_size
(sigma, sigmap[, beam, lattice_reference, eta=0]) Set the beamsize of a particular beam (if the name is given) or all beams (if no name is given) at a reference point called lattice_reference. If lattice_reference is not given the reference point is assumed to be ‘x0’ from the beam in question. This function translates the sizes into twiss parameters. If eta (dispersion) is not given it is assumed to be zero.
Parameters: sigma : list[2]
One sigma for the spatial distribution in the horizontal and vertical directions [sigma_h, sigma_v].
sigmap : list[2]
One sigma for the angular distribution in the horizontal and vertical directions [sigmap_h, sigmap_v].
beam : str
Name of the beam for which these paramters apply.
lattice_reference : list[3]
3D coordinate of there reference point for the known beam dimensions.
eta : float
The dispersion term. Assumed to be zero if not given.
Returns: None
-
oscars.sr.
set_seed
(n) Set the internal random seed
Parameters: n : float
Seed number you wish to use
Returns: None
-
oscars.sr.
set_twiss_parameters
([beam, beta, alpha, gamma, lattice_reference]) Set the twiss parameters for a given beam (if a beam name is given) or for all beams (if no name is given). You should only set 2 of the three [beta, alpha, gamma] as they are not independent.
Parameters: beam : str
Name of the beam for which to set the twiss parameters
beta : list[2]
[horizontal, vertical]
alpha : list[2]
[horizontal, vertical]
gamma : list[2]
[horizontal, vertical]
lattice_reference : list[3]
Location of the values of the twiss parameters. These must be places on the beam axis (though it is not checked).
Returns: None
-
oscars.sr.
version
() Version ID
Returns: version : str
-
oscars.sr.
write_bfield
(oformat[, ofile, bofile, xlim, nx, ylim, ny, zlim, nz, comment, version]) Write the field to ofile in the format described in oformat. You must specify at least one limits and one number of points (e.g. zlim and nz). All output is in SI units.
The formats (oformat) available are: * OSCARS * OSCARS1D * SPECTRA * SRW
For OSCARS1D you also must specify what you want in the output. One spatial dimension must be specified along with at least one field dimension (in any order you like). For examples theses are all valid: ‘OSCARS1D Z Fx Fy Fz’, ‘OSCARS1D Fy Fx Z Fz’.
Parameters: oformat : str
Format of output file
ofile : str
Name of output file
bofile : str
Name of binary output file
xlim : list
[min, max] for x dimension
ylim : list
[min, max] for y dimension
zlim : list
[min, max] for z dimension
comment : str
Comment string to be added to file header. LF and CR are removed.
version : int
Which version of output format (you should use the default unless you have good reason)
Returns: None
-
oscars.sr.
write_efield
(oformat[, ofile, bofile, xlim, nx, ylim, ny, zlim, nz, comment, version]) Write the field to ofile in the format described in oformat. You must specify at least one limits and one number of points (e.g. zlim and nz). All output is in SI units.
The formats (oformat) available are: * OSCARS * OSCARS1D * SPECTRA * SRW
For OSCARS1D you also must specify what you want in the output. One spatial dimension must be specified along with at least one field dimension (in any order you like). For examples theses are all valid: ‘OSCARS1D Z Fx Fy Fz’, ‘OSCARS1D Fy Fx Z Fz’.
Parameters: oformat : str
Format of output file
ofile : str
Name of output file
bofile : str
Name of binary output file
xlim : list
[min, max] for x dimension
ylim : list
[min, max] for y dimension
zlim : list
[min, max] for z dimension
comment : str
Comment string to be added to file header. LF and CR are removed.
version : int
Which version of output format (you should use the default unless you have good reason)
Returns: None
oscars.th¶
OSCARSTH module extension.
-
oscars.th.
add_particle_beam
([type, name, energy_GeV, d0, x0, beam, sigma_energy_GeV, t0, current, weight, rotations, translation, horizontal_direction, beta, alpha, gamma, emittance, eta, lattice_reference, mass, charge]) Add a particle beam to the OSCARS object with a name given by name. There is no limit to the number of different particle beams one can add. They are added with a weight which is by default 1. The weight is used in random sampling when asking for a new particle, for example in oscars.sr.set_new_particle(). If the beam parameter is given you only need to specify name and x0.
- Supported particle types for type are:
- electron
- positron
- muon
- anti-muon
- proton
- anti-proton
- pi+
- pi-
Parameters: type : str
One of the built-in types of particle beams, or ‘custom’. If you use custom you must also specify mass and charge.
name : str
User identified of this beam
energy_GeV : float
Beam energy in [GeV]
d0 : list
Vector [float, float, float] representing the default direction of the beam at the initial position. The normalization of this vector does not matter.
x0 : list
Coordinates of the initial position in [m] [x, y, z]
beam : str
Name of predefined beam
sigma_energy_GeV : float
Beam energy spread in [GeV]
t0 : float
Initial time in [m] at the initial_position. Time here is in terms of ct.
current : float
Beam current in [A]. If this parameter is 0, the current defaults to the single particle charge
weight : float
Weight to give this beam for random sampling when there are multiple beams
rotations : list, optional
3-element list representing rotations around x, y, and z axes: [\(\theta_x, \theta_y, \theta_z\)]
translation : list, optional
3-element list representing a translation in space [x, y, z]
horizontal_direction : list
A list representing the horizontal beam direction. This must be perpendicular to the beam direction. The vertical direction is defined internally as \([direction] \times [horizontal\_direction]\).
beta : list
Values of the horizontal and vertical beta funtion at the point lattice_center [beta_x, beta_y]
alpha : list
Values of the horizontal and vertical alpha funtion at the point lattice_center [alpha_x, alpha_y]
gamma : list
Values of the horizontal and vertical gamma funtion at the point lattice_center [gamma_x, gamma_y]
emittance : list
values of the horizontal and vertical emittance [emittance_x, emittance_y]
eta : list
Values of the horizontal and vertical dispersion at the point lattice_center [eta_x, eta_y]. Currently eta_y is ignored and dispersion is only used in oscars.th calculations
lattice_reference : list
Coordinates of the lattice center [x, y, z] (must be on-axis with respect to the beam)
mass : float
mass of a custom particle. This is only used if type = ‘custom’. Must be non-zero.
charge : float
Charge of a custom particle. This is only used if type = ‘custom’. Must be non-zero.
Returns: None
Currently the predefined beam types are:
- NSLSII
- NSLSII-ShortStraight
- NSLSII-LongStraight
Examples
Add an electron beam with 0.500 [A] current at an initial position of [0, 0, 0] in the Y direction with energy of 3 [GeV]
>>> osr.add_particle_beam(type='electron', name='beam_0', x0=[0, 0, 0], d0=[0, 1, 0], energy_GeV=3, current=0.500)
Add a positron beam with 0.500 [A] current at an initial position of [-2, 0, 0] in the direction given by theta in the X-Y plane with energy of 3 [GeV]
>>> from math import sin, cos >>> theta = 0.25 * osr.pi() >>> osr.add_particle_beam(type='positron', name='beam_0', x0=[-2, 0, 0], d0=[sin(theta), cos(theta), 0], energy_GeV=3, current=0.500)
Add a predefined beam for the NSLSII short straight section
>>> osr.add_particle_beam(beam='NSLSII-ShortStraight', name='beam_0', x0=[-2, 0, 0])
-
oscars.th.
bessel_j
(nu, x) Get the value of the modified bessel function of the second kind K_nu
Parameters: nu : float
x : float
Returns: besselk : float
Value of the bessel function for nu at x
-
oscars.th.
bessel_k
(nu, x) Get the value of the modified bessel function of the second kind K_nu
Parameters: nu : float
x : float
Returns: besselk : float
Value of the bessel function for nu at x
-
oscars.th.
dipole_brightness
() Not Implemented
-
oscars.th.
dipole_critical_energy
(bfield) Get the critical energy for bending magnet in [eV]
Parameters: bfield : float
Magnetic field of dipole in [T]
Returns: energy : float
Dipole critical energy in [eV]
-
oscars.th.
dipole_critical_wavelength
(bfield) Get the critical wavelength for bending magnet in [m]
Parameters: bfield : float
Magnetic field of dipole in [T]
Returns: energy : float
Dipole critical wavelength in [m]
-
oscars.th.
dipole_spectrum
(bfield[, energy_range_eV, energy_points_eV, energy_eV, angle_integrated, angle_range, angle_points, angle, npoints, minimum, ofile, bofile]) Get the spectrum from ideal dipole field. One can calculate the spectrum as a function of photon energy at a given angle or the angular dependence for a particular energy.
Parameters: bfield : float
Magnetic field of the dipole in [T]
energy_range_eV : list
[min, max] photon energy of interest in [eV]
energy_points_eV : list
List of energy points in [eV]
energy_eV : float
Photon energy of interest
angle_integrated : boolean
Set to true if you want the vertical angle integrated
angle_range : list
[min, max] angle of interest in [rad]
angle_points : list
List of angle points in [rad]
angle : float
Angle of interest in [rad]
npoints : int
Number of points for _range requests
ofile : str
Output file name
bofile : str
Binary output file name
Returns: flux : list
A list of flux values for given points [[p0, f0], [p1, f1], ...]. The points can be in energy or angle depending on the input parameters. Output is in [photons / s / mrad^2 / 0.1%bw]
Examples
Calculate the spectrum in a given energy range
>>> oth.dipole_spectrum(bfield=0.4, energy_range_eV=[10, 30000], npoints=1000)
Calculate the spectrum in a given energy range at a vertical angle of 0.0005 [rad]
>>> oth.dipole_spectrum(bfield=0.4, energy_range_eV=[10, 30000], npoints=1000, angle=0.0005)
Calculate the flux in a given angle range at an energy of 2394 [eV]
>>> oth.dipole_spectrum(bfield=0.4, angle_range=[-0.0005, 0.0005], npoints=1000, energy_eV=2394)
-
oscars.th.
print_all
() Print information about what is stored in the th object
Returns: None
-
oscars.th.
print_particle_beams
() Print information about what is stored in the th object
Returns: None
-
oscars.th.
set_particle_beam
([type, name, energy_GeV, d0, x0, beam, sigma_energy_GeV, t0, current, weight, rotations, translation, horizontal_direction, beta, alpha, gamma, emittance, eta, lattice_reference, mass, charge]) Add a particle beam to the OSCARS object with a name given by name. There is no limit to the number of different particle beams one can add. They are added with a weight which is by default 1. The weight is used in random sampling when asking for a new particle, for example in oscars.sr.set_new_particle(). If the beam parameter is given you only need to specify name and x0.
- Supported particle types for type are:
- electron
- positron
- muon
- anti-muon
- proton
- anti-proton
- pi+
- pi-
Parameters: type : str
One of the built-in types of particle beams, or ‘custom’. If you use custom you must also specify mass and charge.
name : str
User identified of this beam
energy_GeV : float
Beam energy in [GeV]
d0 : list
Vector [float, float, float] representing the default direction of the beam at the initial position. The normalization of this vector does not matter.
x0 : list
Coordinates of the initial position in [m] [x, y, z]
beam : str
Name of predefined beam
sigma_energy_GeV : float
Beam energy spread in [GeV]
t0 : float
Initial time in [m] at the initial_position. Time here is in terms of ct.
current : float
Beam current in [A]. If this parameter is 0, the current defaults to the single particle charge
weight : float
Weight to give this beam for random sampling when there are multiple beams
rotations : list, optional
3-element list representing rotations around x, y, and z axes: [\(\theta_x, \theta_y, \theta_z\)]
translation : list, optional
3-element list representing a translation in space [x, y, z]
horizontal_direction : list
A list representing the horizontal beam direction. This must be perpendicular to the beam direction. The vertical direction is defined internally as \([direction] \times [horizontal\_direction]\).
beta : list
Values of the horizontal and vertical beta funtion at the point lattice_center [beta_x, beta_y]
alpha : list
Values of the horizontal and vertical alpha funtion at the point lattice_center [alpha_x, alpha_y]
gamma : list
Values of the horizontal and vertical gamma funtion at the point lattice_center [gamma_x, gamma_y]
emittance : list
values of the horizontal and vertical emittance [emittance_x, emittance_y]
eta : list
Values of the horizontal and vertical dispersion at the point lattice_center [eta_x, eta_y]. Currently eta_y is ignored and dispersion is only used in oscars.th calculations
lattice_reference : list
Coordinates of the lattice center [x, y, z] (must be on-axis with respect to the beam)
mass : float
mass of a custom particle. This is only used if type = ‘custom’. Must be non-zero.
charge : float
Charge of a custom particle. This is only used if type = ‘custom’. Must be non-zero.
Returns: None
Currently the predefined beam types are:
- NSLSII
- NSLSII-ShortStraight
- NSLSII-LongStraight
Examples
Add an electron beam with 0.500 [A] current at an initial position of [0, 0, 0] in the Y direction with energy of 3 [GeV]
>>> osr.add_particle_beam(type='electron', name='beam_0', x0=[0, 0, 0], d0=[0, 1, 0], energy_GeV=3, current=0.500)
Add a positron beam with 0.500 [A] current at an initial position of [-2, 0, 0] in the direction given by theta in the X-Y plane with energy of 3 [GeV]
>>> from math import sin, cos >>> theta = 0.25 * osr.pi() >>> osr.add_particle_beam(type='positron', name='beam_0', x0=[-2, 0, 0], d0=[sin(theta), cos(theta), 0], energy_GeV=3, current=0.500)
Add a predefined beam for the NSLSII short straight section
>>> osr.add_particle_beam(beam='NSLSII-ShortStraight', name='beam_0', x0=[-2, 0, 0])
-
oscars.th.
undulator_K
(bfield, period) Get the undulator K parameter vaule
Parameters: bfield : float
Peak magnetic field in [T]
period : float
Magnetic period in [m]
Returns: K : float
Undulator deflectrion parameter
-
oscars.th.
undulator_bfield
(K, period) Get the undulator BFieldMax parameter value given K and the period
Parameters: K : float
Undulator deflection parameters
period : float
Magnetic period in [m]
Returns: bfield : float
Peak bfield
-
oscars.th.
undulator_brightness
(period, nperiods, harmonic[, bfield_range, K_range, npoints, bfield_points, K_points, minimum, ofile, bofile]) Get the brightness for an ideal undulator given K for a specific harmonic. Should specify either K or bfield, but not both. You must have previously defined a beam, including the beta and emittance values.
Parameters: period : float
Undulator period length [m]
nperiods : int
Number of periods
harmonic : int
Harmonic number of interest
bfield_range : list
[min, max] of bfield range of interest
K_range : list
[min, max] of K range of interest
npoints : int
number of points to use when bfield_range or K_range is specified
bfield_points : list
List of bfield points to calculate brightness at
K_points : list
List of K points to calculate brightness at
minimum : float
Any brightness below the minimum will not be included in the return list
ofile : str
Output file name
bofile : str
Binary output file name
Returns: [energy_eV, brightness]s : list[[float, float], ...]
Photon energy [eV] and brightness [photons/s/0.1%bw/mrad^2/mm] for the given parameters
Examples
Get the harmonic peak spectrum for a undulator with a period of 0.050 [m] having 41 periods and a bfield ranging from 0.05 to 0.8 [T] from the 1st harmonic
>>> oth.undulator_flux_onaxis(period=0.050, nperiods=41, harmonic=1, bfield_range=[0.05, 0.8], npoints=1000)
Get the photon energy and flux for a undulator with a period of 0.050 [m] having 41 periods and K value of 1.8674577 from the 1st harmonic
>>> oth.undulator_flux_onaxis(period=0.050, nperiods=41, harmonic=1, K=1.8674577)
-
oscars.th.
undulator_energy_harmonic
(period, harmonic[, K, bfield]) Get the undulator photon energy at particular harmonic. Specify either K or bfield, but not both.
- period : float
- Undulator period length [m]
- harmonic : int
- Harmonic number of interest
- K : float
- Undulator deflection parameter
- bfield : float
- Magnetic field of interest [T]
Returns: energy : float
Photon energy [eV]
Examples
Find the energy of the 5th harmonic for an undulator of period 0.050 [m] and magnetic field of 0.4 [T]
>>> oth.undulator_energy_harmonic(period=0.050, harmonic=5, bfield=0.4)
Find the energy of the 5th harmonic for an undulator of period 0.050 [m] and K value of 1.8674577
>>> oth.undulator_energy_harmonic(period=0.050, harmonic=5, K=1.8674577)
-
oscars.th.
undulator_flux_onaxis
(period, nperiods, harmonic[, bfield_range, K_range, npoints, bfield_points, K_points, minimum, ofile, bofile]) Get the on-axis flux for an ideal undulator given K for a specific harmonic. Should specify either K or bfield, but not both. You must have previously defined a beam.
Parameters: period : float
Undulator period length [m]
nperiods : int
Number of periods
harmonic : int
Harmonic number of interest
bfield_range : list
[min, max] of bfield range of interest
K_range : list
[min, max] of K range of interest
npoints : int
number of points to use when bfield_range or K_range is specified
bfield_points : list
List of bfield points to calculate flux at
K_points : list
List of K points to calculate flux at
minimum : float
Any flux below the minimum will not be included in the return list
ofile : str
Output file name
bofile : str
Binary output file name
Returns: [energy_eV, flux]s : list[[float, float], ...]
Photon energy [eV] and flux [photons/s/0.1%bw/mrad^2] for the given parameters
Examples
Get the harmonic peak spectrum for a undulator with a period of 0.050 [m] having 41 periods and a bfield ranging from 0.05 to 0.8 [T] from the 1st harmonic
>>> oth.undulator_flux_onaxis(period=0.050, nperiods=41, harmonic=1, bfield_range=[0.05, 0.8], npoints=1000)
Get the photon energy and flux for a undulator with a period of 0.050 [m] having 41 periods and K value of 1.8674577 from the 1st harmonic
>>> oth.undulator_flux_onaxis(period=0.050, nperiods=41, harmonic=1, K=1.8674577)
-
oscars.th.
undulator_period
(bfield, K) Get the undulator Period parameter vaule given this bfield and K
Parameters: K : float
Undulator deflection parameters
bfield : float
Peak magnetic field in [T]
Returns: period : float
Magnetic period
-
oscars.th.
version
() Version ID
Returns: version : str
-
oscars.th.
wiggler_flux
(plane, energy_eV, period, nperiods, width, npoints[, bfield, K, ofile, bofile]) Get the flux for an ideal wiggler according to R. P. Walker XXX XXX. Should specify either K or bfield, but not both. You must have previously defined a beam.
Parameters: period : float
Magnetic period length [m]
nperiods : int
Number of periods
harmonic : int
Harmonic number of interest
bfield_range : list
[min, max] of bfield range of interest
K_range : list
[min, max] of K range of interest
npoints : int
number of points to use when bfield_range or K_range is specified
bfield_points : list
List of bfield points to calculate flux at
K_points : list
List of K points to calculate flux at
minimum : float
Any flux below the minimum will not be included in the return list
ofile : str
Output file name
bofile : str
Binary output file name
Returns: [energy_eV, flux]s : list[[float, float], ...]
Photon energy [eV] and flux [photons/s/0.1%bw/mrad^2] for the given parameters
Examples
Get the harmonic peak spectrum for a undulator with a period of 0.050 [m] having 41 periods and a bfield ranging from 0.05 to 0.8 [T] from the 1st harmonic
>>> oth.undulator_flux_onaxis(period=0.050, nperiods=41, harmonic=1, bfield_range=[0.05, 0.8], npoints=1000)
Get the photon energy and flux for a undulator with a period of 0.050 [m] having 41 periods and K value of 1.8674577 from the 1st harmonic
>>> oth.undulator_flux_onaxis(period=0.050, nperiods=41, harmonic=1, K=1.8674577)
-
oscars.th.
wiggler_spectrum
(bfield, period, length[, energy_range_eV, energy_points_eV, energy_eV, angle_integrated, angle_range, angle_points, angle, npoints, minimum, ofile, bofile]) Get the spectrum from ideal dipole field. One can calculate the spectrum as a function of photon energy at a given angle or the angular dependence for a particular energy.
Parameters: bfield : float
Magnetic field of the dipole in [T]
period : float
Length of magnetic period in [m]
length : float
Length of device in [m]
energy_range_eV : list
[min, max] photon energy of interest in [eV]
energy_points_eV : list
List of energy points in [eV]
energy_eV : float
Photon energy of interest
angle_integrated : boolean
Set to true if you want the vertical angle integrated
angle_range : list
[min, max] angle of interest in [rad]
angle_points : list
List of angle points in [rad]
angle : float
Angle of interest in [rad]
npoints : int
Number of points for _range requests
ofile : str
Output file name
bofile : str
Binary output file name
Returns: flux : list
A list of flux values for given points [[p0, f0], [p1, f1], ...]. The points can be in energy or angle depending on the input parameters. Output is in [photons / s / mrad^2 / 0.1%bw]
Examples
Calculate the spectrum in a given energy range
>>> oth.wiggler_spectrum(bfield=0.4, period=0.100, length=7.0, energy_range_eV=[10, 30000], npoints=1000)
oscars.util¶
-
oscars.util.
add_power_density
(p)[source]¶ Add power densities from list of power densities and return new summed power density
-
oscars.util.
read_file_list
(ifile, idir=None)[source]¶ read a list of parameters and filenames from a file.
The file format should be two columns separated by whitespace. In the first column is the parameter value as a floating point number and in the second column the file name. Whitespace in the name is unforgivably wrong, but technically allowed.
Parameters: ifile : str
Full path to file containing the list of parameters and filenames
idir : str
Path to directory where the files are contained if not in the ‘ifile’ directory
Returns: file_list : list
List of [parameter, filename]s
-
oscars.util.
read_file_list2
(ifile, gap=None, phase_mode=None, phase=None, idir=None)[source]¶ read a list of parameters and filenames from a file.
The file format should be four columns separated by whitespace. The columns are gap value, phase mode, phase value, file name. Whitespace in the phase mode is not allowed.
- One can specify the following:
- gap + phase_mode phase + phase_mode
Parameters: ifile : str
Full path to file containing the list of parameters and filenames
gap : float
Gap value of interest (must be exact match, not interpolated)
phase_mode : str
Phase mode of interest
phase : float
Phase value
idir : str
Path to directory where the files are contained if not in the ‘ifile’ directory
Returns: file_list : list
List of [gap, phase_mode, phase, filename]s
oscars.plots_mpl¶
-
oscars.plots_mpl.
plot_bfield
(osr, mymin=-1, mymax=1, ylim=None, show=True, ofile='', axis='Z', npoints=20000, between_two_points=None, ret=False)[source]¶ Plot the magnetic field as a function of Z
-
oscars.plots_mpl.
plot_efield
(osr, mymin=-1, mymax=1, ylim=None, show=True, ofile='', axis='Z', npoints=20000, between_two_points=None, ret=False)[source]¶ Plot the electric field as a function of Z
-
oscars.plots_mpl.
plot_electric_field_vs_time
(efield, show=True, ofile='', ret=False)[source]¶ Plot the electric field as a function of time
-
oscars.plots_mpl.
plot_flux
(V, title='Flux [$\\gamma / mm^2 / 0.1\\%bw / s$]', xlabel='X1 Axis [$m$]', ylabel='X2 Axis [$m$]', clim=None, show=True, ofile='', figsize=None, ylim=None, xlim=None, colorbar=True, ret=False, nticks_cb=None)[source]¶ Plot a 2D histogram with equal spacing
-
oscars.plots_mpl.
plot_flux_spectrum
(F, S, energy=None, title='Flux [$\\gamma / mm^2 / 0.1\\%bw / s]$', xlabel='X1 Axis [$m$]', ylabel='X2 Axis [$m$]', show=True, ofile='', figsize=[10, 3], ylim=None, xlim=None, colorbar=True, ret=False)[source]¶ Plot a 2D histogram with equal spacing
-
oscars.plots_mpl.
plot_power_density
(V, title=None, xlabel='X1 Axis [$m$]', ylabel='X2 Axis [$m$]', show=True, ofile='', figsize=None, ret=False, x1=None, x2=None)[source]¶ Plot a 2D histogram with equal spacing
-
oscars.plots_mpl.
plot_power_density_1d
(V, title='Power Density [$W / mm^2$]', xlabel='[$m$]', ylabel='[$W / mm^2$]', xlim=None, ylim=None, show=True, ofile='', figsize=None, ret=False)[source]¶ Plot a 1D power density
-
oscars.plots_mpl.
plot_power_density_2d1d
(V, x1=None, x2=None, title=None, xlabel='[$m$]', ylabel='[$W / mm^2$]', xlim=None, ylim=None, show=True, ofile='', figsize=None, ret=False)[source]¶ Plot a 2D histogram with equal spacing
-
oscars.plots_mpl.
plot_spectrum
(S, log=False, show=True, ofile='', title='Spectrum', xlabel='Energy [eV]', ylabel='$\\gamma / mm^2 / 0.1\\%bw / s$', figsize=None, ylim=None, xlim=None, transparent=True, ret=False, xticks=None, **kwargs)[source]¶ Plot the spectrum
-
oscars.plots_mpl.
plot_trajectory_position
(trajectory, show=True, ofile='', axis='Z', figsize=[18, 4.5], ret=False)[source]¶ Plot the trajectory position. You can optionally change the axis, output to a file, show or not, and return or not
Parameters: - trajectory (list) – Particle trajectory
- show (bool) – to show the plot or not
- ofile (str) – output file name
- axis (str) – which axis to plot trajectory on
- figsize (list) – dimension of the figure
- ret (bool) – to return the plot or not
-
oscars.plots_mpl.
plot_trajectory_velocity
(trajectory, show=True, ofile='', figsize=[18, 4.5], ret=False)[source]¶ Plot the trajectory velocity. You can optionally change the axis, output to a file, show or not, and return or not
Parameters: - trajectory (list) – Particle trajectory
- show (bool) – to show the plot or not
- ofile (str) – output file name
- axis (str) – which axis to plot trajectory on
- figsize (list) – dimension of the figure
- ret (bool) – to return the plot or not
-
oscars.plots_mpl.
plot_undulator_brightness
(oth, period, nperiods, harmonics, minimum=0, bfield=None, K=None, show=True, ret=False, title='Brightness', figsize=None, ofile=None)[source]¶ Plot the brightness of an undulator.
oscars.plots3d_mpl¶
-
oscars.plots3d_mpl.
plot_bfield2D
(srs, xlim=[-0.001, 0.001], zlim=[-1, 1], nx=20, nz=50)[source]¶ plot the bfield in 3D vector form
-
oscars.plots3d_mpl.
plot_bfield3D
(srs, xlim=[-0.02, 0.02], ylim=[-0.02, 0.02], zlim=[-0.2, 0.02], nx=10, ny=10, nz=10)[source]¶ plot the bfield in 3D vector form
-
oscars.plots3d_mpl.
plot_surface
(surface, xlim=None, ylim=None, zlim=None, **kwargs)[source]¶ plot a parametric surface in 3d
-
oscars.plots3d_mpl.
power_density_3d
(srs, surface, normal=1, rotations=[0, 0, 0], translation=[0, 0, 0], nparticles=0, gpu=0, nthreads=0, title='Power Density [$W / mm^2$]', xlim=None, ylim=None, zlim=None, colorbar=True, figsize=None, alpha=0.4, ofile=None, show=True, view_init=None, axis=None, transparent=True, xticks=None, yticks=None, zticks=None, bbox_inches='tight', max_level=24, quantity='power density')[source]¶ calculate power density for and plot a parametric surface in 3d
oscars.parametric_surfaces¶
-
class
oscars.parametric_surfaces.
PSCone
(r0=0, r1=1, L=1, ustart=None, ustop=None, nu=10, vstart=None, vstop=None, nv=10)[source]¶ A Parametric surface - cylinder with no top or bottom
Methods
normal
(u, v)Return a unit normal in 3D at this u and v position position
(u, v)Return the position in 3D at this u and v
-
class
oscars.parametric_surfaces.
PSCurvedRectangle
(L=1, W=1, R=0, ustart=-0.5, ustop=0.5, nu=10, vstart=-0.5, vstop=0.5, nv=10)[source]¶ A Parametric surface - curved rectangle
Methods
normal
(u, v)Return a unit normal in 3D at this u and v position position
(u, v)Return the position in 3D at this u and v
-
class
oscars.parametric_surfaces.
PSCylinder
(R=1, L=1, ustart=None, ustop=None, nu=10, vstart=None, vstop=None, nv=10)[source]¶ A Parametric surface - cylinder with no top or bottom
Methods
normal
(u, v)Return a unit normal in 3D at this u and v position position
(u, v)Return the position in 3D at this u and v
-
class
oscars.parametric_surfaces.
PSDisk
(r0=1, r1=1, nu=10, vstart=0, vstop=6.283185307179586, nv=10)[source]¶ A Parametric surface - Disk
Methods
normal
(u, v)Return a unit normal in 3D at this u and v position position
(u, v)Return the position in 3D at this u and v
-
class
oscars.parametric_surfaces.
PSRectangle
(L=1, W=1, ustart=-0.5, ustop=0.5, nu=10, vstart=-0.5, vstop=0.5, nv=10)[source]¶ A Parametric surface - rectangle
Methods
normal
(u, v)Return a unit normal in 3D at this u and v position position
(u, v)Return the position in 3D at this u and v
-
class
oscars.parametric_surfaces.
PSSphere
(R=1, ustart=0, ustop=6.283185307179586, nu=10, vstart=-1.5707963267948966, vstop=1.5707963267948966, nv=10)[source]¶ A Parametric surface - sphere
Methods
normal
(u, v)Return a unit normal in 3D at this u and v position position
(u, v)Return the position in 3D at this u and v
-
class
oscars.parametric_surfaces.
PSTorus
(R=1, r=1, ustart=0, ustop=6.283185307179586, nu=10, vstart=0, vstop=6.283185307179586, nv=10)[source]¶ A Parametric surface - torus
Methods
normal
(u, v)Return a unit normal in 3D at this u and v position position
(u, v)Return the position in 3D at this u and v