ocelot.optics.wave

wave optics

Module Contents

Classes

RadiationField

3d or 2d coherent radiation distribution, *.fld variable is the same as Genesis dfl structure

WaistScanResults

TransferFunction

data container for Fourier Optics transfer functions

StokesParameters

HeightProfile

1d surface of mirror

Spectrogram

spectrogram of the pulse

WignerDistribution

calculated Wigner distribution (spectrogram) of the pulse

Functions

bin_stokes(S, bin_size)

needs fix for 3d!!!

calc_stokes_out(out1, out2, pol='rl', on_axis=True)

calc_stokes_dfl(dfl1, dfl2, basis='xy', mode=(0, 0))

mode: (average_longitudinally, sum_transversely)

calc_stokes(E1, E2, s=None, basis='xy')

average_stokes_l(S, sc_range=None)

sum_stokes_tr(S)

generate_dfl(*args, **kwargs)

generate_gaussian_dfl(xlamds=1e-09, shape=(51, 51, 100), dgrid=(0.001, 0.001, 5e-05), power_rms=(0.0001, 0.0001, 5e-06), power_center=(0, 0, None), power_angle=(0, 0), power_waistpos=(0, 0), wavelength=None, zsep=None, freq_chirp=0, en_pulse=None, power=1000000.0, **kwargs)

generates RadiationField object

imitate_sase_dfl(xlamds, rho=0.0002, seed=None, **kwargs)

imitation of SASE radiation in 3D

calc_phase_delay_poly(coeff, w, w0)

Calculate the phase delay with given coefficients

screen2dfl(screen, polarization='x')

Function converts synchrotron radiation from ocelot.rad.screen.Screen to ocelot.optics.wave.RadiationField.

dfl_disperse(dfl, coeff, E_ph0=None, return_result=False)

The function adds a phase shift to the fld object in the frequency domain

dfl_ap(*args, **kwargs)

dfl_ap_rect(dfl, ap_x=np.inf, ap_y=np.inf)

model rectangular aperture to the radaition in either domain

dfl_ap_circ(dfl, r=np.inf, center=(0, 0))

apply circular aperture to the radaition in either domain

dfl_prop(dfl, z, fine=1, debug=1)

LEGACY, WILL BE DEPRECATED, SEE METHOD

dfl_waistscan(dfl, z_pos, projection=0, **kwargs)

propagates the RadaitionField object dfl

dfl_interp(dfl, interpN=(1, 1), interpL=(1, 1), newN=(None, None), newL=(None, None), method='cubic', return_result=1, **kwargs)

2d interpolation of the coherent radiation distribution

dfl_shift_z(dfl, s, set_zeros=1, return_result=1)

shift the radiation within the window in time domain

dfl_pad_z(dfl, padn, return_result=1)

dfl_cut_z(dfl, z=[-np.inf, np.inf], debug=1)

dfl_fft_z(dfl, method='mp', nthread=multiprocessing.cpu_count(), debug=1)

LEGACY, WILL BE DEPRECATED, SEE METHOD

dfl_fft_xy(dfl, method='mp', nthread=multiprocessing.cpu_count(), debug=1)

LEGACY, WILL BE DEPRECATED, SEE METHOD

dfl_trf(dfl, trf, mode, dump_proj=False)

Multiplication of radiation field by given transfer function (transmission or ferlection, given by mode)

trf_mult(trf_list, embed_list=True)

multiply transfer functions

calc_phase_delay(coeff, w, w0)

expression for the phase – coeff[0] + coeff[1]*(w - w0)/1! + coeff[2]*(w - w0)**2/2! + coeff[3]*(w - w0)**3/3!

dfl_chirp_freq(dfl, coeff, E_ph0=None, return_result=False)

The function adds a phase shift to a fld object. The expression for the phase see in the calc_phase_delay function

generate_1d_profile(hrms, length=0.1, points_number=1000, wavevector_cutoff=0, psd=None, seed=None)

Function for generating HeightProfile of highly polished mirror surface

dfl_reflect_surface(dfl, angle, hrms=None, height_profile=None, axis='x', seed=None, return_height_profile=False)

Function models the reflection of ocelot.optics.wave.RadiationField from the mirror surface considering effects

trf_mult_mix(trf_list, mode_out='ref')

multiply transfer functions in a mixed way:

save_trf(trf, attr, flePath)

calc_wigner(field, method='mp', nthread=multiprocessing.cpu_count(), debug=1)

calculation of the Wigner distribution

wigner_pad(wig, pad)

pads WignerDistribution with zeros in time domain

wigner_out(out, z=inf, method='mp', pad=1, debug=1, on_axis=1)

returns WignerDistribution from GenesisOutput at z

wigner_dfl(dfl, method='mp', pad=1, **kwargs)

returns on-axis WignerDistribution from dfl file

wigner_stat(out_stat, stage=None, z=inf, method='mp', debug=1, pad=1, **kwargs)

returns averaged WignerDistribution from GenStatOutput at stage at z

wigner_smear(wig, sigma_s)

Convolves wigner distribution with gaussian window function to obtain spectrogram, see https://arxiv.org/pdf/1811.11446.pdf

calc_ph_sp_dens(spec, freq_ev, n_photons, spec_squared=1)

calculates number of photons per electronvolt

imitate_1d_sase_like(td_scale, td_env, fd_scale, fd_env, td_phase=None, fd_phase=None, phen0=None, en_pulse=None, fit_scale='td', n_events=1, **kwargs)

Models FEL pulse(s) based on Gaussian statistics

imitate_1d_sase(spec_center=500, spec_res=0.01, spec_width=2.5, spec_range=(None, None), pulse_length=6, en_pulse=0.001, flattop=0, n_events=1, spec_extend=5, **kwargs)

Models FEL pulse(s) based on Gaussian statistics

dfldomain_check(domains, both_req=False)

ocelot.optics.wave._logger
ocelot.optics.wave.nthread
ocelot.optics.wave.fftw_avail = True
ocelot.optics.wave.__author__ = Svitozar Serkez, Andrei Trebushinin, Mykola Veremchuk
class ocelot.optics.wave.RadiationField(shape=0, 0, 0)

3d or 2d coherent radiation distribution, *.fld variable is the same as Genesis dfl structure

fileName(self)
copy_param(self, dfl1, version=1)
__getitem__(self, i)
__setitem__(self, i, fld)
shape(self)

returns the shape of fld attribute

domains(self)

returns domains of the radiation field

Lz(self)

full longitudinal mesh size

Ly(self)

full transverse vertical mesh size

Lx(self)

full transverse horizontal mesh size

Nz(self)

number of points in z

Ny(self)

number of points in y

Nx(self)

number of points in x

intensity(self)

3d intensity, abs(fld)**2

int_z(self)

intensity projection on z power [W] or spectral density [arb.units]

ang_z_onaxis(self)

on-axis phase

int_y(self)

intensity projection on y

int_x(self)

intensity projection on x

int_xy(self)
int_zx(self)
int_zy(self)
E(self)

energy in the pulse [J]

scale_kx(self)
scale_ky(self)
scale_kz(self)
scale_x(self)
scale_y(self)
scale_z(self)
ph_sp_dens(self)
curve_wavefront(self, r=np.inf, plane='xy', domain_z=None)

introduction of the additional wavefront curvature with radius r

r can be scalar or vector with self.Nz() points r>0 -> converging wavefront

plane is the plane in which wavefront is curved:

‘x’ - horizontal focusing ‘y’ - vertical focusing ‘xy’ - focusing in both planes

domain_z is the domain in which wavefront curvature is introduced

‘f’ - frequency ‘t’ - time None - original domain (default)

to_domain(self, domains='ts', **kwargs)

tranfers radiation to specified domains *domains is a string with one or two letters:

(“t” or “f”) and (“s” or “k”)

where

‘t’ (time); ‘f’ (frequency); ‘s’ (space); ‘k’ (inverse space);

e.g.

‘t’; ‘f’; ‘s’; ‘k’; ‘ts’; ‘fs’; ‘tk’; ‘fk’

order does not matter

**kwargs are passed down to self.fft_z and self.fft_xy

fft_z(self, method='mp', nthread=multiprocessing.cpu_count(), **kwargs)
fft_xy(self, method='mp', nthread=multiprocessing.cpu_count(), **kwargs)
prop(self, z, fine=1, return_result=0, return_orig_domains=1, **kwargs)

Angular-spectrum propagation for fieldfile

can handle wide spectrum
(every slice in freq.domain is propagated

according to its frequency)

no kx**2+ky**2<<k0**2 limitation

dfl is the RadiationField() object z is the propagation distance in [m] fine=1 is a flag for ~2x faster propagation.

no Fourier transform to frequency domain is done assumes no angular dispersion (true for plain FEL radiation) assumes narrow spectrum at center of xlamds (true for plain FEL radiation)

return_result does not modify self, but returns result

z>0 -> forward direction

prop_m(self, z, m=1, fine=1, return_result=0, return_orig_domains=1, **kwargs)

Angular-spectrum propagation for fieldfile

can handle wide spectrum
(every slice in freq.domain is propagated

according to its frequency)

no kx**2+ky**2<<k0**2 limitation

dfl is the RadiationField() object z is the propagation distance in [m] m is the output mesh size in terms of input mesh size (m = L_out/L_inp) which can be a number m or a pair of number m = [m_x, m_y] fine==0 is a flag for ~2x faster propagation.

no Fourier transform to frequency domain is done assumes no angular dispersion (true for plain FEL radiation) assumes narrow spectrum at center of xlamds (true for plain FEL radiation)

z>0 -> forward direction

mut_coh_func(self, norm=1, jit=1)

calculates mutual coherence function consider downsampling the field first

coh(self, jit=0)

calculates degree of transverse coherence consider downsampling the field first

tilt(self, angle=0, plane='x', return_orig_domains=True)

deflects the radaition in given direction by given angle by introducing transverse phase chirp

disperse(self, disp=0, E_ph0=None, plane='x', return_orig_domains=True)

introducing angular dispersion in given plane by deflecting the radaition by given angle depending on its frequency disp is the dispertion coefficient [rad/eV] E_ph0 is the photon energy in [eV] direction of which would not be changed (principal ray)

class ocelot.optics.wave.WaistScanResults
fileName(self)
class ocelot.optics.wave.TransferFunction

data container for Fourier Optics transfer functions

ev(self)
__mul__(self, f)
class ocelot.optics.wave.StokesParameters
__getitem__(self, i)
P_pol(self)
P_pol_l(self)
deg_pol(self)
deg_pol_l(self)
deg_pol_c(self)
chi(self)
psi(self)
slice_2d(self, loc, plane='z')
slice_2d_idx(self, idx, plane='z')
proj(self, plane='x', mode='sum')
class ocelot.optics.wave.HeightProfile

1d surface of mirror

hrms(self)
set_hrms(self, rms)
psd(self)
save(self, path)
load(self, path, return_obj=False)
ocelot.optics.wave.bin_stokes(S, bin_size)

needs fix for 3d!!!

ocelot.optics.wave.calc_stokes_out(out1, out2, pol='rl', on_axis=True)
ocelot.optics.wave.calc_stokes_dfl(dfl1, dfl2, basis='xy', mode=0, 0)

mode: (average_longitudinally, sum_transversely)

ocelot.optics.wave.calc_stokes(E1, E2, s=None, basis='xy')
ocelot.optics.wave.average_stokes_l(S, sc_range=None)
ocelot.optics.wave.sum_stokes_tr(S)
class ocelot.optics.wave.Spectrogram

spectrogram of the pulse (always positive!)

__init(self)
class ocelot.optics.wave.WignerDistribution

calculated Wigner distribution (spectrogram) of the pulse in time/frequency domain as space/wavelength

property freq_lamd(self)
power(self)
spectrum(self)
energy(self)
fileName(self)
eval(self, method='mp')
inst_freq(self)
group_delay(self)
inst_bandwidth(self)
ocelot.optics.wave.generate_dfl(*args, **kwargs)
ocelot.optics.wave.generate_gaussian_dfl(xlamds=1e-09, shape=51, 51, 100, dgrid=0.001, 0.001, 5e-05, power_rms=0.0001, 0.0001, 5e-06, power_center=0, 0, None, power_angle=0, 0, power_waistpos=0, 0, wavelength=None, zsep=None, freq_chirp=0, en_pulse=None, power=1000000.0, **kwargs)

generates RadiationField object narrow-bandwidth, paraxial approximations

xlamds [m] - central wavelength shape (x,y,z) - shape of field matrix (reversed) to dfl.fld dgrid (x,y,z) [m] - size of field matrix power_rms (x,y,z) [m] - rms size of the radiation distribution (gaussian) power_center (x,y,z) [m] - position of the radiation distribution power_angle (x,y) [rad] - angle of further radiation propagation power_waistpos (x,y) [m] downstrean location of the waist of the beam wavelength [m] - central frequency of the radiation, if different from xlamds zsep (integer) - distance between slices in z as zsep*xlamds freq_chirp dw/dt=[1/fs**2] - requency chirp of the beam around power_center[2] en_pulse, power = total energy or max power of the pulse, use only one

ocelot.optics.wave.imitate_sase_dfl(xlamds, rho=0.0002, seed=None, **kwargs)

imitation of SASE radiation in 3D

xlamds - wavelength of the substracted fast-varying component rho - half of the expected FEL bandwidth **kwargs identical to generate_dfl()

returns RadiationField object

ocelot.optics.wave.calc_phase_delay_poly(coeff, w, w0)

Calculate the phase delay with given coefficients coeff — coefficients in phase delay expression: w — photon frequencies w0 — photon frequency at which zero group delay is introduced

The expression for the phase: delta_phi = coeff[0] +

coeff[1]*(w - w0)**1/1!/(1e15)**1 + coeff[2]*(w - w0)**2/2!/(1e15)**2 + coeff[3]*(w - w0)**3/3!/(1e15)**3 + …

… coeff[n]*(w - w0)**n/n!/(1e15)**n

coeff is an array-like object with coeff[0] — phase measured in [rad] coeff[1] — group delay measured in [fs ^ 1] coeff[2] — group delay dispersion (GDD) measured in [fs ^ 2] coeff[3] — third-order dispersion (TOD) measured in [fs ^ 3] … coeff[n] — nth-order dispersion measured in [fs ^ n]

@author Andrei Trebushinin

ocelot.optics.wave.screen2dfl(screen, polarization='x')

Function converts synchrotron radiation from ocelot.rad.screen.Screen to ocelot.optics.wave.RadiationField. New ocelot.optics.wave.RadiationField object will be generated without changing ocelot.rad.screen.Screen object.

Parameters
  • screen – ocelot.rad.screen.Screen object, electric field of which will be used to generate RadiationField

  • polarization – polarization for conversion to RadiationField (‘x’ or ‘y’)

Returns

ocelot.optics.wave.RadiationField in domains = ‘fs’

ocelot.optics.wave.dfl_disperse(dfl, coeff, E_ph0=None, return_result=False)

The function adds a phase shift to the fld object in the frequency domain

dfl — is a RadiationField object coeff — coefficients in phase delay expression:

The expression for the phase: delta_phi = coeff[0] +

coeff[1]*(w - w0)**1/1!/(1e15)**1 + coeff[2]*(w - w0)**2/2!/(1e15)**2 + coeff[3]*(w - w0)**3/3!/(1e15)**3 + …

… coeff[n]*(w - w0)**n/n!/(1e15)**n

coeff is an array-like object with coeff[0] — phase measured in [rad] coeff[1] — group delay measured in [fs ^ 1] coeff[2] — group delay dispersion (GDD) measured in [fs ^ 2]

e.g. coeff[2]: (dt[fs] / dw[1/fs]) = 1e-6 * 1e15**2 * hr_eV_s / speed_of_light * (ds [um] / dE [eV])

coeff[3] — third-order dispersion (TOD) measured in [fs ^ 3] … coeff[n] — nth-order dispersion measured in [fs ^ n]

can be
  • a list (coeff[0], coeff[1], … coeff[n])

  • or a number, then coeff[2] is expected

E_ph0 — photon energy at which zero group delay is introduced
can be:
  • None - when the carrier freq will be extracted from dfl object

  • center - when the carrier freq will be calculated as averaging over time

  • from dfl.int_z() and dfl.scale_z()

  • E_ph0 - a number in eV

return_result — a flag that is responsible for returning the dfl object
if return_result is True:

create a copy and return the modified copy without affecting the original object

else:

change the original dfl object

@author Andrei Trebushinin

ocelot.optics.wave.dfl_ap(*args, **kwargs)
ocelot.optics.wave.dfl_ap_rect(dfl, ap_x=np.inf, ap_y=np.inf)

model rectangular aperture to the radaition in either domain

ocelot.optics.wave.dfl_ap_circ(dfl, r=np.inf, center=0, 0)

apply circular aperture to the radaition in either domain

ocelot.optics.wave.dfl_prop(dfl, z, fine=1, debug=1)

LEGACY, WILL BE DEPRECATED, SEE METHOD

Fourier propagator for fieldfile

can handle wide spectrum
(every slice in freq.domain is propagated

according to its frequency)

no kx**2+ky**2<<k0**2 limitation

dfl is the RadiationField() object z is the propagation distance in [m] fine==0 is a flag for ~2x faster propagation.

no Fourier transform to frequency domain is done assumes no angular dispersion (true for plain FEL radiation) assumes narrow spectrum at center of xlamds (true for plain FEL radiation)

returns RadiationField() object

z>0 ==> forward

ocelot.optics.wave.dfl_waistscan(dfl, z_pos, projection=0, **kwargs)

propagates the RadaitionField object dfl through the sequence of positions z_pos and calculates transverse distribution parameters such as peak photon density and sizes in both dimensions

if projection==1, then size of projection is calculated

otherwise - size across the central line passing through the mesh center

ocelot.optics.wave.dfl_interp(dfl, interpN=1, 1, interpL=1, 1, newN=None, None, newL=None, None, method='cubic', return_result=1, **kwargs)

2d interpolation of the coherent radiation distribution interpN and interpL define the desired interpolation coefficients for transverse point __density__ and transverse mesh __size__ correspondingly newN and newL define the final desire number of points and size of the mesh when newN and newL are not None interpN and interpL values are ignored coordinate convention is (x,y)

ocelot.optics.wave.dfl_shift_z(dfl, s, set_zeros=1, return_result=1)

shift the radiation within the window in time domain dfl - initial RadiationField object s - longitudinal offset value in meters set_zeros - to set the values outside the time window to zeros

ocelot.optics.wave.dfl_pad_z(dfl, padn, return_result=1)
ocelot.optics.wave.dfl_cut_z(dfl, z=[- np.inf, np.inf], debug=1)
ocelot.optics.wave.dfl_fft_z(dfl, method='mp', nthread=multiprocessing.cpu_count(), debug=1)

LEGACY, WILL BE DEPRECATED, SEE METHOD

ocelot.optics.wave.dfl_fft_xy(dfl, method='mp', nthread=multiprocessing.cpu_count(), debug=1)

LEGACY, WILL BE DEPRECATED, SEE METHOD

ocelot.optics.wave.dfl_trf(dfl, trf, mode, dump_proj=False)

Multiplication of radiation field by given transfer function (transmission or ferlection, given by mode) dfl is RadiationField() object (will be mutated) trf is TransferFunction() object mode is either ‘tr’ for transmission

or ‘ref’ for reflection

ocelot.optics.wave.trf_mult(trf_list, embed_list=True)

multiply transfer functions trf_list is a list of transfer functions embed_list == True will write the list of input transfer functions into the output transfer function as an trf.trf_list instance

returns TransferFunction() object

ocelot.optics.wave.calc_phase_delay(coeff, w, w0)

expression for the phase – coeff[0] + coeff[1]*(w - w0)/1! + coeff[2]*(w - w0)**2/2! + coeff[3]*(w - w0)**3/3! coeff is a list with coeff[0] =: measured in [rad] — phase coeff[1] =: measured in [fm s ^ 1] — group delay coeff[2] =: measured in [fm s ^ 2] — group delay dispersion (GDD) coeff[3] =: measured in [fm s ^ 3] — third-order dispersion (TOD) …

ocelot.optics.wave.dfl_chirp_freq(dfl, coeff, E_ph0=None, return_result=False)

The function adds a phase shift to a fld object. The expression for the phase see in the calc_phase_delay function dfl — is a fld object coeff — coefficients in phase (see in the calc_phase_delay function) E_ph0 — energy with respect to which the phase shift is calculated return_result — a flag that is responsible for returning the modified dfl object if it is True or

change the dfl function parameter if it is False

ocelot.optics.wave.generate_1d_profile(hrms, length=0.1, points_number=1000, wavevector_cutoff=0, psd=None, seed=None)

Function for generating HeightProfile of highly polished mirror surface

Parameters
  • hrms – [meters] height errors root mean square

  • length – [meters] length of the surface

  • points_number – number of points (pixels) at the surface

  • wavevector_cutoff – [1/meters] point on k axis for cut off small wavevectors (large wave lengths) in the PSD (with default value 0 effects on nothing)

  • psd – [meters^3] 1d array; power spectral density of surface (if not specified, will be generated) (if specified, must have shape = (points_number // 2 + 1, ), otherwise it will be cut to appropriate shape)

  • seed – seed for np.random.seed() to allow reproducibility

Returns

HeightProfile object

ocelot.optics.wave.dfl_reflect_surface(dfl, angle, hrms=None, height_profile=None, axis='x', seed=None, return_height_profile=False)

Function models the reflection of ocelot.optics.wave.RadiationField from the mirror surface considering effects of mirror surface height errors. The method based on phase modulation. The input RadiationField object is modified

Parameters
  • dfl – RadiationField object from ocelot.optics.wave

  • angle – [radians] angle of incidence with respect to the surface

  • hrms – [meters] height root mean square of reflecting surface

  • height_profile – HeightProfile object of the reflecting surface (if not specified, will be generated using hrms)

  • axis – direction along which reflection takes place

  • seed – seed for np.random.seed() to allow reproducibility

  • return_height_profile – boolean type variable; if it equals True the function will return height_profile

ocelot.optics.wave.trf_mult_mix(trf_list, mode_out='ref')

multiply transfer functions in a mixed way: trf_list is list of tulpes, like [(trf1,’ref’),(trf2,’tr’)], here ‘ref’ and ‘tr’ mean that reflectivity trom transter function trf1 is multiplied by transmissivity of transfer function trf2 mode_out is a string ‘ref’ or ‘tr’ that specifies into thich instance to write the multiplied output embed_list == True will write the list of input transfer functions into the output transfer function as an trf.trf_list instance

returns TransferFunction() object

ocelot.optics.wave.save_trf(trf, attr, flePath)
ocelot.optics.wave.calc_wigner(field, method='mp', nthread=multiprocessing.cpu_count(), debug=1)

calculation of the Wigner distribution input should be an amplitude and phase of the radiation as list of complex numbers with length N output is a real value of wigner distribution

ocelot.optics.wave.wigner_pad(wig, pad)

pads WignerDistribution with zeros in time domain

ocelot.optics.wave.wigner_out(out, z=inf, method='mp', pad=1, debug=1, on_axis=1)

returns WignerDistribution from GenesisOutput at z

ocelot.optics.wave.wigner_dfl(dfl, method='mp', pad=1, **kwargs)

returns on-axis WignerDistribution from dfl file

ocelot.optics.wave.wigner_stat(out_stat, stage=None, z=inf, method='mp', debug=1, pad=1, **kwargs)

returns averaged WignerDistribution from GenStatOutput at stage at z

ocelot.optics.wave.wigner_smear(wig, sigma_s)

Convolves wigner distribution with gaussian window function to obtain spectrogram, see https://arxiv.org/pdf/1811.11446.pdf

Parameters
  • wig – ocelot.optics.wave.WignerDistribution object which will be convolved with generated window func

  • sigma_s – [meters] rms size of the s=-ct the gaussian window func

Returns

convolved ocelot.optics.wave.WignerDistribution object with gaussian window function

ocelot.optics.wave.calc_ph_sp_dens(spec, freq_ev, n_photons, spec_squared=1)

calculates number of photons per electronvolt

ocelot.optics.wave.imitate_1d_sase_like(td_scale, td_env, fd_scale, fd_env, td_phase=None, fd_phase=None, phen0=None, en_pulse=None, fit_scale='td', n_events=1, **kwargs)

Models FEL pulse(s) based on Gaussian statistics td_scale - scale of the pulse on time domain [m] td_env - expected pulse envelope in time domain [W] fd_scale - scale of the pulse in frequency domain [eV] fd_env - expected pulse envelope in frequency domain [a.u.] td_phase - additional phase chirp to be added in time domain fd_phase - additional phase chirp to be added in frequency domain phen0 - sampling wavelength expressed in photon energy [eV] en_pulse - expected average energy of the pulses [J] fit_scale - defines the scale in which outputs should be returned:

‘td’ - time domain scale td_scale is used for the outputs, frequency domain phase and envelope will be re-interpolated ‘fd’ - frequency domain scale fd_scale is used for the outputs, time domain phase and envelope will be re-interpolated

n_events - number of spectra to be generated

returns tuple of 4 arguments: (ph_en, fd, s, td) fd_scale - colunm of photon energies in eV fd - matrix of radiation in frequency domain with shape, normalized such that np.sum(abs(fd)**2) is photon spectral density, i.e: np.sum(abs(fd)**2)*fd_scale = N_photons td - matrix of radiation in time domain, normalized such that abs(td)**2 = radiation_power in [w]

ocelot.optics.wave.imitate_1d_sase(spec_center=500, spec_res=0.01, spec_width=2.5, spec_range=None, None, pulse_length=6, en_pulse=0.001, flattop=0, n_events=1, spec_extend=5, **kwargs)

Models FEL pulse(s) based on Gaussian statistics spec_center - central photon energy in eV spec_res - spectral resolution in eV spec_width - width of spectrum in eV (fwhm of E**2) spec_range = (E1, E2) - energy range of the spectrum. If not defined, spec_range = (spec_center - spec_width * spec_extend, spec_center + spec_width * spec_extend) pulse_length - longitudinal size of the pulse in um (fwhm of E**2) en_pulse - expected average energy of the pulses in Joules flattop - if true, flat-top pulse in time domain is generated with length ‘pulse_length’ in um n_events - number of spectra to be generated

return tuple of 4 arguments: (s, td, ph_en, fd) ph_en - colunm of photon energies in eV with size (spec_range[2]-spec_range[1])/spec_res fd - matrix of radiation in frequency domain with shape ((spec_range[2]-spec_range[1])/spec_res, n_events), normalized such that np.sum(abs(fd)**2) is photon spectral density, i.e: np.sum(abs(fd)**2)*spec_res = N_photons s - colunm of longitudinal positions along the pulse in yime domain in um td - matrix of radiation in time domain with shape ((spec_range[2]-spec_range[1])/spec_res, n_events), normalized such that abs(td)**2 = radiation_power

ocelot.optics.wave.dfldomain_check(domains, both_req=False)