Member Functions

Primer for the Examples

All examples here assume that you have the OSCARS SR module in your path and have an OSCARS SR object called osr, perhaps as in the following example

# Import the OSCARS SR module
import oscars.sr

# Creat an OSCARS SR object
osr = oscars.sr.sr()

oscars.sr

oscars.sr.pi()

Get the value of pi: 3.14159265358979323846

Returns:float
oscars.sr.qe()

Get the value of elementary charge: +1.602176462e-19 [C]

Returns:float
oscars.sr.rand()

Random float between 0 and 1

Returns:float
oscars.sr.norm()

Random float from the normal distribution

Returns:float
oscars.sr.set_seed(seed)

Set the internal random seed

Parameters:seed (int) – Random seed
oscars.sr.set_gpu_global(gpu)

If set to 1, OSCARS will try to use the GPU for all calculations

Parameters:gpu (int) – 1 or 0
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:integer
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 multi-threading.

oscars.sr.get_ctstart()

Gets the start time for particle propogation. The time is in units of [m].

Returns:float
oscars.sr.get_ctstop()

Gets the stop time for particle propogation. The time is in units of [m].

Returns:float
oscars.sr.set_ctstartstop(time_start, time_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:
  • time_start (float) – Start time
  • time_stop (float) – Stop time
Returns:

None

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.get_npoints_trajectory()

Gets the number of points to be used in the trajectory calculation

Returns:int - Number of points
oscars.sr.add_bfield_file(ifile, iformat[, rotation, translation, scaling])

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 scaling 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\)

scaling 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
  • iformat (str) – Format of the input file
  • rotation (list[3]) – [\(\theta_x, \theta_y, \theta_z\)]
  • translation – Translation in space [x, y, z]
  • scaling – Scaling of input parameters
Returns:

None

Example:
# Add a magnetic field from a file where the columns are in the order Z, Bx, By, Bz where Z is in [m] and Bx, By, Bz are in [T].
osr.add_bfield_file(ifile='file.txt', iformat='OSCARS1D Z Bx By Bz')
oscars.sr.add_bfield_function(function)

Adds a magnetic field in the form of a user defined python function. The input for this function must be (x, y, z, t).

Parameters:function (func) – Python function
Returns:None
Example:
# Create a function in python and use it as a magnetic 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])

Add a gaussian magnetic 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 the following: \(\\exp(-(x - x_0)^2 / \\sigma_x^2)\)

Parameters:
  • bfield (list[3]) – A list representing the peak field [Bx, By, Bz]
  • center (list[3]) – A list representing the coordinates of the center point [X0, Y0, Z0]
  • sigma (list[3]) – A list representing the sigma in each coordinate [S0, S0, S0]
Returns:

None

Example:

Add a magnetic field of 0.5 [T] in the X-direction centered at X=0, Y=0, Z=0 [in meters], with a sigma of 0.1 [m] in the Z-direction

# Will add a magnetic field of 1 Tesla in the Y-direction centered
# at X=0, Y=0, Z=0 [in meters], with a sigma of 0.1 in the Z-direction
osr.add_bfield_gaussian(bfield=[1, 0, 0], sigma=[0, 0, 0.10])
Example:Add a magnetic field of 1 Tesla in the Y-direction centered at X=1, Y=1, Z=1 [in meters], with a sigma of 0.05 in the Z-direction
# Will add a magnetic field of 1 Tesla in the Y-direction centered
# at X=1, Y=1, Z=1 [in meters], 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_uniform(bfield[, width, rotations, translation])

Add a uniform magnetic 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[3]) – A list representing the magnetic field [Bx, By, Bz]
  • width (list[3]) – A list representing the spetial extent of the field [Wx, Wy, Wz]
  • rotations (list[3]) – A list representing the rotations of this object about the X, Y amd Z axis (in that order)
  • translation (list[3]) – A list representing the translation of the center of this object
Returns:

None

Example:

Add a magnetic field of 0.0001 [T] in the X-direction for all space

# Will add a magnetic field of 0.0001 [T] in the X-direction over all space
osr.add_bfield_uniform(bfield=[0.0001, 0, 0])
Example:Add a magnetic 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].
# Will add a magnetic field of 1 Tesla in the Y-direction centered
# at X=0, Y=0, Z=0.75 [in meters], extending in Z a width of 1.5 [m]
osr.add_bfield_uniform(bfield=[0, 1, 0], width=[0, 0, 1.5], translation=[0, 0, 0.75])
oscars.sr.add_bfield_undulator(bfield, period, nperiods[, phase, rotations, translation])

Adds an ideal sinusoidal undulator field with a given maximum b-field amplitude, period, and number of periods. Optionally one can specify the phase offset (in [rad]), rotations and translation.

Parameters:
  • bfield (list[3]) – A list representing the peak field [Bx, By, Bz] in [T]
  • period (float) – Length of one period
  • nperiods (int) – Number of periods
  • phase (float) – Phase offset in [rad].
  • rotations (list[3]) – A list representing the rotations of this object about the X, Y amd Z axis (in that order)
  • translation (list[3]) – A list representing the translation of the center of this object
Returns:

None

Example:

Add the magnetic field for an undulator of length 3 [m] with a period of 0.050 [m] (60 periods) having a peak magnetic field in the Y direction of 1.1 [T]

# Will add magnetic field 3 [m] undulaotor having 1.1 [T] peak field
# in the Y-direction with a period of 0.050 [m]
osr.add_bfield_undulator(bfield=[0, 1.1, 0], period=0.050, nperiods=60)
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[3]) – A 3D list representing a point in space [x, y, z]
Returns:[float, float, float] - list representing the B-field [Bx, By, Bz]
oscars.sr.clear_bfields()

Remove all of the existing magnetic fields.

Returns:None
oscars.sr.add_efield_file(ifile, iformat[, rotation, translation, scaling])

Add an 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 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 scaling list to scale the input (which must be in SI units of [m] for distances/positions and [T] for electric field values.

The rotation is performed first in the order: \(\theta_x, \theta_y, \theta_z\)

scaling 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
  • iformat (str) – Format of the input file
  • rotation (list[3]) – [\(\theta_x, \theta_y, \theta_z\)]
  • translation – Translation in space [x, y, z]
  • scaling – Scaling of input parameters
Returns:

None

Example:
# Add electric field from a file where the columns are in the order Z, Bx, By, Bz where Z is in [m] and Bx, By, Bz are in [T].
osr.add_efield_file(ifile='file.txt', iformat='OSCARS1D Z Bx By Bz')
oscars.sr.add_efield_function(function)

Adds electric field in the form of a user defined python function. The input for this function must be (x, y, z, t).

Parameters:function (func) – Python function
Returns:None
Example:
# Create a function in python and use it as a electric 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 electric 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 the following: \(\\exp(-(x - x_0)^2 / \\sigma_x^2)\)

Parameters:
  • efield (list[3]) – A list representing the peak field [Bx, By, Bz]
  • center (list[3]) – A list representing the coordinates of the center point [X0, Y0, Z0]
  • sigma (list[3]) – A list representing the sigma in each coordinate [S0, S0, S0]
Returns:

None

Example:

Add a electric field of 0.5 [T] in the X-direction centered at X=0, Y=0, Z=0 [in meters], with a sigma of 0.1 [m] in the Z-direction

# Will add a electric field of 1 Tesla in the Y-direction centered
# at X=0, Y=0, Z=0 [in meters], with a sigma of 0.1 in the Z-direction
osr.add_efield_gaussian(efield=[1, 0, 0], sigma=[0, 0, 0.10])
Example:Add a electric field of 1 Tesla in the Y-direction centered at X=1, Y=1, Z=1 [in meters], with a sigma of 0.05 in the Z-direction
# Will add a electric field of 1 Tesla in the Y-direction centered
# at X=1, Y=1, Z=1 [in meters], 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_uniform(efield[, width, rotations, translation])

Add a uniform electric 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[3]) – A list representing the electric field [Bx, By, Bz]
  • width (list[3]) – A list representing the spetial extent of the field [Wx, Wy, Wz]
  • rotations (list[3]) – A list representing the rotations of this object about the X, Y amd Z axis (in that order)
  • translation (list[3]) – A list representing the translation of the center of this object
Returns:

None

Example:

Add a electric field of 0.0001 [T] in the X-direction for all space

# Will add a electric field of 0.0001 [T] in the X-direction over all space
osr.add_efield_uniform(efield=[0.0001, 0, 0])
Example:Add a electric 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].
# Will add a electric field of 1 Tesla in the Y-direction centered
# at X=0, Y=0, Z=0.75 [in meters], extending in Z a width of 1.5 [m]
osr.add_efield_uniform(efield=[0, 1, 0], width=[0, 0, 1.5], translation=[0, 0, 0.75])
oscars.sr.add_efield_undulator(efield, period, nperiods[, phase, rotations, translation])

Adds an ideal sinusoidal undulator field with a given maximum b-field amplitude, period, and number of periods. Optionally one can specify the phase offset (in [rad]), rotations and translation.

Parameters:
  • efield (list[3]) – A list representing the peak field [Bx, By, Bz] in [T]
  • period (float) – Length of one period
  • nperiods (int) – Number of periods
  • phase (float) – Phase offset in [rad].
  • rotations (list[3]) – A list representing the rotations of this object about the X, Y amd Z axis (in that order)
  • translation (list[3]) – A list representing the translation of the center of this object
Returns:

None

Example:

Add the electric field for an undulator of length 3 [m] with a period of 0.050 [m] (60 periods) having a peak electric field in the Y direction of 1.1 [T]

# Will add electric field 3 [m] undulaotor having 1.1 [T] peak field
# in the Y-direction with a period of 0.050 [m]
osr.add_efield_undulator(efield=[0, 1.1, 0], period=0.050, nperiods=60)
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[3]) – A 3D list representing a point in space [x, y, z]
Returns:[float, float, float] - list representing the B-field [Bx, By, Bz]
oscars.sr.clear_efields()

Remove all of the existing electric fields.

Returns:None
oscars.sr.write_bfield(ofile, oformat[, xlim, nx, ylim, ny, zlim, nz, comment])

Write the magnetic 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 Bx By Bz’, ‘OSCARS1D By Bx Z Bz’.

Parameters:
  • ofile (str) – Name of output file
  • oformat (str) – format of output file
  • xlim (list[min, max]) – min and max for x dimension
  • ylim (list[min, max]) – min and max for y dimension
  • zlim (list[min, max]) – min and max for z dimension
  • comment (str) – comment string to be added to file header. LF and CR are removed.
Returns:

None

oscars.sr.write_efield(ofile, oformat[, xlim, nx, ylim, ny, zlim, nz, comment])

Same as write_bfield, but writes the electric field.

oscars.sr.set_particle_beam(type, name, energy_GeV, d0, x0[, sigma_energy_GeV, t0, current, weight, rotations, translation, horizontal_direction, beta, emittance, lattice_center, mass, charge])

This function is the same as oscars.sr.add_particle_beam(), but it clears all particle beams before the ‘add’.

oscars.sr.add_particle_beam(type, name, energy_GeV, x0, d0[, sigma_energy_GeV, t0, current, weight, rotations, translation, horizontal_direction, beta, emittance, 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().

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 ([float, float, float]) – Vector representing the default direction of the beam at the initial position. The normalization of this vector does not matter.
  • x0 ([float, float, float]) – Coordinates of the initial position in [m]
  • sigma_energy_GeV (float) – Beam energy 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[3]) – A list representing the rotations of this beam about the X, Y amd Z axis (in that order)
  • translation (list[3]) – A list representing the translation of the x0 of this object
  • horizontal_direction (list[3]) – 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 – values of the horizontal and vertical beta funtion at the point lattice_center [beta_x, beta_y]
  • emittance ([float, float]) – values of the horizontal and vertical emittance [emittance_x, emittance_y]
  • lattice_reference (list[3]) – Coordinates of the lattice center (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

Example:

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]

# Will add an electron beam called beam_0 starting at [0, 0, 0] headed in the
# direction of +Y [0, 1, 0], with an energy of 3 [GeV] at time T0 = 0 [m] having a
# current of 0.500 [A] and a weight 1
osr.add_particle_beam(type='electron', name='beam_0', x0=[0, 0, 0], d0=[0, 1, 0], energy_GeV=3, current=0.500)
Example: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]
# Create a positron beam in the X-Y plane at an angle theta
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)
oscars.sr.clear_particle_beams()

Remove all of the existing particle beams

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.get_particle_x0()

Get the initial position for the current particle in [m]

Returns:[float, float, float]
oscars.sr.get_particle_beta0()

Get the initial \(\vec \beta\) for the current particle

Returns:[float, float, float]
oscars.sr.get_particle_e0()

Get the initial energy for the current particle in [GeV]

Returns:float
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.

Returns:A list of points of the form [[[x, y, z], [Beta_x, Beta_y, Beta_z]], ...]
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]], ...]

Returns:A list of points of the form [[[x, y, z], [Beta_x, Beta_y, Beta_z]], ...]
oscars.sr.calculate_spectrum(obs[, npoints, energy_range_eV, points_eV, polarization, horizontal_direction, propogation_direction, angle, nparticles, ofile, nthreads, gpu])

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 ([float, float, float]) – Point where you wish to calculate the spectrum
  • npoints (int) – Number of points to calculate in the given energy range
  • energy_range_eV ([float, float]) – energy range in eV as a list of length 2
  • points_eV ([float, ..]) – A list of points to calculate the flux at
  • 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 ([x, y, z]) – The direction you consider to be horizontal. Should be perpendicular to the photon beam propogation direction
  • vertical_direction – Same as horizontal_direction but the vertical direction
  • 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
  • ofile (str) – File name you would like to output this data to
  • nthreads (int) – Number of threads to use
  • gpu (int) – Use the gpu or not
Returns:

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], ...]

Example:

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.

# 0.4 [T] 2 [m] long dipole centered at Z=0 [m]
osr.add_bfield_uniform(bfield=[0, 0.4, 0], width=[0, 0, 2])

# NSLS2 electron beam
osr.add_particle_beam(type='electron', name='beam_0', x0=[0, 0, 0], direction=[0, 0, 1], t0=0, energy_GeV=3, current=0.500)

# Set start and stop time for calculation
osr.set_ctstartstop(-0.5, 0.5)


# Set number of points for trajectory calculation (hight for dipole)
osr.set_npoints_trajectory(20000)

# Calculate spectrum zt Z=30 [m] in the energy range of 100 to
# 1000 [eV] at 900 equally spaced 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:float - Total power in [W]
oscars.sr.calculate_power_density_rectangle(npoints[, plane, width, x0x1x2, rotations, translation, ofile, normal, dim, nthreads, gpu])

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, 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 ([float, float]) – Width if rectangle in X1 and X2
  • x0x1x2 ([[float, float, float], [float, float, float], [float, float, float]]) – List of three points [[x0, y0, z0], [x1, y1, z1], [x2, y2, z2]]
  • rotations (list[3]) – A list representing the rotations of this beam about the X, Y amd Z axis (in that order)
  • translation (list[3]) – A list representing the translation of the x0 of this object
  • ofile (str) – 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.
  • 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.
  • nthreads (int) – Number of threads to use
  • gpu (int) – Use the gpu or not
Returns:

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], ...]

Example:

Calculate the power density within a simple rectangle 1 [cm] x 1 [cm], 30 [m] downstream

# Calculate power density in rectangle in the XY plane 30 [m] downstream
# from a photon beam in the +Z direction
power_density = osr.calculate_power_density(plane='XY', width=[0.01, 0.01], npoints=[51, 51], translation=[0, 0, 30])
Example:Calculate the power density within a simple rectangle 1 [cm] x 1 [cm], 30 [m] downstream defined using the three-point method
# Calculate power density in rectangle in the XY plane 30 [m] downstream
# from a photon beam in the +Z direction
power_density = osr.calculate_power_density(x0x1x2=[[-0.005, -0.005, 30], [+0.005, -0.005, 30], [-0.005, +0.005, 30]], npoints=[51, 51])
Example:Calculate the power density on a flat surface close to and parallel to the beam direction
# Calculate power density in rectangle on the upper flat portion of a beampipe in the middle
# of an undulator.  Here the 'normal' is reversed to get the correct sign
osr.add_undulator(bfield=[0, 1, 0], period=0.050, nperiods=40)
power_density = osr.calculate_power_density(plane='XZ', width=[0.008, 2], npoints=[51, 101], normal=-1)
oscars.sr.calculate_power_density(points[, normal, rotations, translation, ofile, nthreads, gpu])

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[ list[list[float, float, float], list[float, float, float]], ..]) – 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).
  • 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[3]) – A list representing the rotations of this beam about the X, Y amd Z axis (in that order)
  • translation (list[3]) – A list representing the translation of the x0 of this object
  • ofile (str) – Output file name
  • nthreads (int) – Number of threads to use
  • gpu (int) – Use the gpu or not
Returns:

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], ...]

oscars.sr.calculate_flux_rectangle(energy_eV, npoints[, plane, width, x0x1x2, rotations, translation, polarization, horizontal_direction, propogation_direction, angle, ofile, normal, dim, nthreads, gpu])

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) – The photon energy you are interested in
  • npoints ([int, 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 ([float, float]) – Width if rectangle in X1 and X2
  • x0x1x2 ([[float, float, float], [float, float, float], [float, float, float]]) – List of three points [[x0, y0, z0], [x1, y1, z1], [x2, y2, z2]]
  • rotations (list[3]) – A list representing the rotations of this beam about the X, Y amd Z axis (in that order)
  • translation (list[3]) – A list representing the translation of the x0 of this object
  • 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 ([x, y, z]) – The direction you consider to be horizontal. Should be perpendicular to the photon beam propogation direction
  • vertical_direction – Same as horizontal_direction but the vertical direction
  • 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
  • ofile (str) – 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.
  • 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.
  • nthreads (int) – Number of threads to use
  • gpu (int) – Use the gpu or not
Returns:

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], ...]

oscars.sr.average_spectra([ifiles, bifiles, ofile, bofile])

Average spectra from different files and output to specified file. The input files must have the same format. Binary operation will be supported in a future release.

You must specify an input and an output

Parameters:
  • ifiles (list) – The input file names
  • ofile (str) – The output file name
  • bifiles (list) – The binary input file names
  • bofile (str) – The binary output file name
Returns:

None

oscars.sr.average_flux([ifiles, bifiles, ofile, bofile, dim])

Average flux from different files and output to specified file. The input files must have the same format. Binary operation will be supported in a future release.

You must specify an input and an output

Parameters:
  • ifiles (list) – The input file names
  • ofile (str) – The output file name
  • bifiles (list) – The binary input file names
  • bofile (str) – The binary output file name
  • dim (integer) – in 2 or 3 dimensions (default is 2)
Returns:

None

oscars.sr.average_power_density([ifiles, bifiles, ofile, bofile, dim])

Average power densities from different files and output to specified file. The input files must have the same format. Binary operation will be supported in a future release.

You must specify an input and an output

Parameters:
  • ifiles (list) – The input file names
  • ofile (str) – The output file name
  • bifiles (list) – The binary input file names
  • bofile (str) – The binary output file name
  • dim (integer) – in 2 or 3 dimensions (default is 2)
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.get_spectrum()

Get the current spectrum stored in memory

Returns:list
oscars.sr.add_to_flux(flux[, weight])

Add a flux to the current flux with a given weight.

Parameters:
  • spectrum (list) – a list of [[x, y, z], flux] pairs (flux format)
  • weight (float) – Weight for this flux
Returns:

None

oscars.sr.get_flux()

Get the current flux map stored in memory

Returns:list of [[x, y, z], flux]
oscars.sr.add_to_power_density(power_density[, weight])

Add a power density to the current power density with a given weight.

Parameters:
  • spectrum (list) – a list of [[x, y, z], power_density] pairs (power density format)
  • weight (float) – Weight for this power density
Returns:

None

oscars.sr.get_power_density()

Get the current power_density map stored in memory

Returns:list of [[x, y, z], power_density]
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 ([float, float, float]) – Point where you wish to calculate the electric field
  • ofile (str) – Output file name
Returns:

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.plots_mpl

oscars.plots_mpl.plot_bfield(osc, 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(osc, 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$]', show=True, ofile='', figsize=None, ylim=None, xlim=None, colorbar=True, ret=False)[source]

Plot a 2D histogram with equal spacing

oscars.plots_mpl.plot_power_density(V, title='Power Density [$W / mm^2$]', xlabel='X1 Axis [$m$]', ylabel='X2 Axis [$m$]', 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', figsize=None, ylim=None, xlim=None, transparent=True, ret=False, **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.total_power(pd)[source]

Calculate the total power in a uniform grid.

This will not work for a non-uniform grid. Different NX and NY are ok.

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, ret=False, title='Power Density', xlim=None, ylim=None, zlim=None, colorbar=True, figsize=None, alpha=0.4, ofile=None, show=True, view_init=None, axis=None, transparent=True)[source]

calculate power density for and plot a parametric surface in 3d

oscars.parametric_surfaces

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

normal(u, v)[source]

Return a unit normal in 3D at this u and v position

position(u, v)[source]

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

normal(u, v)[source]

Return a unit normal in 3D at this u and v position

position(u, v)[source]

Return the position in 3D at this u and v

class oscars.parametric_surfaces.PSSphere(R=1, ustart=None, ustop=None, nu=10, vstart=None, vstop=None, nv=10)[source]

A Parametric surface - sphere

normal(u, v)[source]

Return a unit normal in 3D at this u and v position

position(u, v)[source]

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

normal(u, v)[source]

Return a unit normal in 3D at this u and v position

position(u, v)[source]

Return the position in 3D at this u and v