qdyn.pulse module

Module containing the Pulse class and functions for initializing pulse shapes.

Summary

Classes:

Pulse

Numerical real or complex control pulse

Functions:

CRAB_carrier

Construct a “carrier” based on the CRAB formula

blackman

Return a Blackman function between t_start and t_stop, see http://en.wikipedia.org/wiki/Window_function#Blackman_windows

carrier

Create the “carrier” of the pulse as a weighted superposition of cosines at different frequencies.

gaussian

Return a Gaussian shape with peak amplitude 1.0

pulse_tgrid

Return a pulse time grid suitable for an equidistant time grid of the states between t0 and T with nt intervals.

tgrid_from_config

Extract the time grid from the given config file

Reference

class qdyn.pulse.Pulse(tgrid, amplitude=None, time_unit=None, ampl_unit=None, freq_unit=None, config_attribs=None)[source]

Bases: object

Numerical real or complex control pulse

Parameters
  • tgrid (numpy.ndarray(float)) – Time grid values

  • amplitude (numpy.ndarray(float), numpy.ndarray(complex)) – Amplitude values. If not given, amplitude will be zero

  • time_unit (str) – Unit of values in tgrid. Will be ignored when reading from file.

  • ampl_unit (str) – Unit of values in amplitude. Will be ignored when reading from file.

  • freq_unit (str) – Unit of frequencies when calculating spectrum. If not given, an appropriate unit based on time_unit will be chosen, if possible (or a TypeError will be raised.

tgrid

time points at which the pulse values are defined, from t0 + dt/2 to T - dt/2.

Type

numpy.ndarray(float)

amplitude

array of real or complex pulse values.

Type

numpy.ndarray(float), numpy.ndarray(complex)

time_unit

Unit of values in tgrid

Type

str

ampl_unit

Unit of values in amplitude

Type

str

freq_unit

Unit to use for frequency when calculating the spectrum

Type

str

preamble

List of lines that are written before the header when writing the pulse to file. Each line should start with ‘# ‘

Type

list

postamble

List of lines that are written after all data lines. Each line should start with ‘# ‘

Type

list

config_attribs

Additional config data, for when generating a QDYN config file section describing the pulse (e.g. {‘oct_shape’: ‘flattop’, ‘t_rise’: ‘10_ns’})

Type

dict

Class Attributes:
unit_convert (QDYN.units.UnitConvert): converter to be used for any

unit conversion within any methods

Example

>>> tgrid = pulse_tgrid(10, 100)
>>> amplitude = 100 * gaussian(tgrid, 50, 10)
>>> p = Pulse(tgrid=tgrid, amplitude=amplitude, time_unit='ns',
...           ampl_unit='MHz')
>>> p.write('pulse.dat')
>>> p2 = Pulse.read('pulse.dat')
>>> from os import unlink; unlink('pulse.dat')

Notes

It is important to remember that the QDYN Fortran library considers pulses to be defined on the intervals of the propagation time grid (i.e. for a time grid with n time steps of dt, the pulse will have n-1 points, defined at points shifted by dt/2)

The pulse_tgrid and tgrid_from_config routine may be used to obtain the proper pulse time grid from the propagation time grid:

>>> import numpy as np
>>> p = Pulse(tgrid=pulse_tgrid(10, 100), ampl_unit='MHz',
...           time_unit='ns')
>>> len(p.tgrid)
99
>>> print(str(p.dt))
0.10101_ns
>>> p.t0
0
>>> print("%.5f" % p.tgrid[0])
0.05051
>>> print(str(p.T))
10_ns
>>> print("%.5f" % p.tgrid[-1])
9.94949

The type of the amplitude (not whether there is a non-zero imaginary part) decide whether the pulse is considered real or complex. Complex pulses are not allowed to couple to Hermitian operators, and in an optimization, both the real and imaginary part of the pulse are modified.

unit_convert = <qdyn.units.UnitConvert object>
freq_units = {'au': 'au', 'dimensionless': 'dimensionless', 'fs': 'eV', 'iu': 'iu', 'microsec': 'MHz', 'ns': 'GHz', 'ps': 'cminv', 'unitless': 'unitless'}
__eq__(other)[source]

Compare two pulses, within a precision of 1e-12

copy()[source]

Return a copy of the pulse

classmethod read(filename, time_unit=None, ampl_unit=None, freq_unit=None, ignore_header=False)[source]

Read a pulse from file, in the format generated by the QDYN write_pulse routine.

Parameters
  • filename (str) – Path and name of file from which to read the pulse

  • time_unit (str or None) – The unit of the time grid

  • ampl_unit (str or None) – The unit of the pulse amplitude.

  • freq_unit (str or None) – Intended value for the freq_unit attribute. If None, a freq_unit will be chosen automatically, if possible (or a TypeError will be raised)

  • ignore_header (bool) – If True, the file header will be ignored.

Note

By default, the file is assumed to contain a header that identifies the columns and their units, as a comment line immediately preceding the data. If time_unit or ampl_unit are None, and ignore_header is False, the respective unites are extracted from the header line. If time_unit or ampl_unit are not None, the respective values will be converted from the unit specified in the file header. If ignore_header is True, both time_unit and ampl_unit must be given. This can be used to read in pulses that were not generated by the QDYN write_pulse routine. Note that if ignore_header is True, all comment lines preceding the data will be included in the preamble attribute.

The write method allows to restore exactly the original pulse file.

classmethod from_func(tgrid, func, time_unit=None, ampl_unit=None, freq_unit=None, config_attribs=None)[source]

Instantiate a pulse from an amplitude function func.

All other parameters are passed on to __init__

property dt

Time grid step, as instance of UnitFloat

property t0

Time at which the pulse begins (dt/2 before the first point in the pulse), as instance of UnitFloat

property states_tgrid

Time grid values for the states propagated under the numerical pulse values, as numpy array in units of time_unit.

The returned time grid has one point more than tgrid, and extends from t0 to T (inclusive).

property w_max

Maximum frequency that can be represented with the current sampling rate.

property dw

Step width in the spectrum (i.e. the spectral resolution) based on the current pulse duration, as an instance of UnitFloat.

property T

Time at which the pulse ends (dt/2 after the last point in the pulse), as an instance of UnitFloat.

property is_complex

Is the pulse amplitude of complex type?

as_func(interpolation='linear', allow_args=False)[source]

Callable that evaluates the pulse for a given time value.

Possible values for interpolation are ‘linear’ and ‘piecewise’.

The resulting function takes an argument t that must be a float in the range [t0, T] and in units of time_unit). It returns the (interpolated) pulse amplitude as a float, in units of ampl_unit

If allow_args is True, the resulting function takes a second argument args that is ignored. This is for compatibility with qutip, see http://qutip.org/docs/latest/guide/dynamics/dynamics-time.html.

convert(time_unit=None, ampl_unit=None, freq_unit=None)[source]

Convert the pulse data to different units

get_timegrid_point(t, move='left')[source]

Return the next point to the left (or right) of the given t which is on the pulse time grid

property fluence

Fluence (integrated pulse energy) for the pulse

\[\int_{-\infty}^{\infty} \vert|E(t)\vert^2 dt\]
property oct_iter

OCT iteration number from the pulse preamble, if available. If not available, 0

spectrum(freq_unit=None, mode='complex', sort=False)[source]

Calculate the spectrum of the pulse

Parameters
  • freq_unit (str) – Desired unit of the freq output array. Can Hz (GHz, Mhz, etc) to obtain frequencies, or any energy unit, using the correspondence f = E/h. If not given, defaults to the freq_unit attribute

  • mode (str) – Wanted mode for spectrum output array. Possible values are ‘complex’, ‘abs’, ‘real’, ‘imag’

  • sort (bool) – Sort the output freq array (and the output spectrum array) so that frequecies are ordered from -w_max .. 0 .. w_max, instead of the direct output from the FFT. This is good for plotting, but does not allow to do an inverse Fourier transform afterwards

Returns

Frequency values associated with the amplitude values in spectrum, i.e. the x-axis of the spectrogram. The values are in the unit freq_unit. Real (mode in [‘abs’, ‘real’, ‘imag’]) or complex (mode=’complex’) amplitude of each frequency component.

Return type

numpy.ndarray(float), numpy.ndarray(complex)

Notes

If sort=False and mode=’complex’, the original pulse values can be obtained by simply calling np.fft.ifft

The spectrum is not normalized (Scipy follows the convention of doing the normalization on the backward transform). You might want to normalized by 1/n for plotting.

fftfreq(freq_unit=None)[source]

Return the FFT frequencies associated with the pulse. Cf. numpy.fft.fftfreq

Parameters

freq_unit (str) – Desired unit of the output array. If not given, defaults to the freq_unit attribute

Returns

Frequency values associated with the pulse time grid. The first half of the freq array contains the positive frequencies, the second half the negative frequencies

Return type

numpy.ndarray(float)

derivative()[source]

Calculate the derivative of the current pulse and return it as a new pulse. Note that the derivative is in units of ampl_unit/time_unit, but will be marked as ‘unitless’.

phase(unwrap=False, s=None, derivative=False, freq_unit=None)[source]

Return the pulse’s complex phase, or derivative of the phase

Parameters
  • unwrap (bool) – If False, report the phase in [-pi:pi]. If True, the phase may take any real value, avoiding the discontinuous jumps introduced by limiting the phase to a 2 pi interval.

  • s (float or None) – smoothing parameter, see scipy.interpolate.UnivariateSpline. If None, no smoothing is performed.

  • derivative (bool) – If False, return the (smoothed) phase directly. If True, return the derivative of the (smoothed) phase.

  • freq_unit (str or None) – If derivative is True, the unit in which the derivative should be calculated. If None, self.freq_unit is used.

Note

When calculating the derivative, some smoothing is generally required. By specifying a smoothing parameter s, the phase is smoothed through univeriate splines before calculating the derivative.

When calculating the phase directly (instead of the derivative), smoothing should only be used when also unwrapping the phase.

write(filename, mode=None)[source]

Write a pulse to file, in the same format as the QDYN write_pulse routine

Parameters
  • filename (str) – Name of file to which to write the pulse

  • mode (str) – Mode in which to write files. Possible values are ‘abs’, ‘real’, or ‘complex’. The former two result in a two-column file, the latter in a three-column file. If not given, ‘real’ or ‘complex’ is used, depending on the type of amplitude

write_oct_spectral_filter(filename, filter_func, freq_unit=None)[source]

Evaluate a spectral filter function and write the result to the file with a given filename, in a format such that the file may be used for the oct_spectral_filter field of a pulse in a QDYN config file. The file will have two columns: The pulse frequencies (see fftfreq method), and the value of the filter function in the range [0, 1]

Parameters
  • filename (str) – Filename of the output file

  • filter_func (callable) – A function that takes a frequency values (in units of freq_unit) and returns a filter value in the range [0, 1]

  • freq_unit (str) – Unit of frequencies that filter_func assumes. If not given, defaults to the freq_unit attribute.

Note

The filter_func function may return any values that numpy considers equivalent to floats in the range [0, 1]. This includes boolean values, where True is equivalent to 1.0 and False is equivalent to 0.0

apply_spectral_filter(filter_func, freq_unit=None)[source]

Apply a spectral filter function to the pulse (in place)

Parameters
  • filter_func (callable) – A function that takes a frequency values (in units of freq_unit) and returns a filter value in the range [0, 1]

  • freq_unit (str) – Unit of frequencies that filter_func assumes. If not given, defaults to the freq_unit attribute.

apply_smoothing(**kwargs)[source]

Smooth the pulse amplitude (in place) through univariate splining. All keyword arguments are passed directly to scipy.interpolate.UnivariateSpline. This especially includes the smoothing parameter s.

resample(upsample=None, downsample=None, num=None, window=None)[source]

Resample the pulse, either by giving an upsample ratio, a downsample ration, or a number of sampling points

Parameters
  • upsample (int) – Factor by which to increase the number of samples. Afterwards, those points extending beyond the original end point of the pulse are discarded.

  • downsample (int) – For downsample=n, keep only every n’th point of the original pulse. This may cause the resampled pulse to end earlier than the original pulse

  • num (int) – Resample with num sampling points. This may case the end point of the resampled pulse to change

  • window (list, numpy.ndarray, callable, str, float, or tuple) – Specifies the window applied to the signal in the Fourier domain. See sympy.signal.resample.

Notes

Exactly one of upsample, downsample, or num must be given.

Upsampling will maintain the pulse start and end point (as returned by the T and t0 properties), up to some rounding errors. Downsampling, or using an arbitrary number will change the end point of the pulse in general.

render_pulse(ax, label='pulse')[source]

Render the pulse amplitude on the given axes.

render_phase(ax, label='phase')[source]

Render the complex phase of the pulse on the given axes.

render_spectrum(ax, zoom=True, wmin=None, wmax=None, spec_scale=None, spec_max=None, freq_unit=None, mark_freqs=None, mark_freq_points=None, label='spectrum')[source]

Render spectrum onto the given axis, see plot for arguments

plot(fig=None, show_pulse=True, show_spectrum=True, zoom=True, wmin=None, wmax=None, spec_scale=None, spec_max=None, freq_unit=None, mark_freqs=None, mark_freq_points=None, **figargs)[source]

Generate a plot of the pulse on a given figure

Parameters
  • fig (matplotlib.figure.Figure) – The figure onto which to plot. If not given, create a new figure from matplotlib.pyplot.figure

  • show_pulse (bool) – Include a plot of the pulse amplitude? If the pulse has a vanishing imaginary part, the plot will show the real part of the amplitude, otherwise, there will be one plot for the absolute value of the amplitude and one showing the complex phase in units of pi

  • show_spectrum (bool) – Include a plot of the spectrum?

  • zoom (bool) – If True, only show the part of the spectrum that has amplitude of at least 0.1% of the maximum peak in the spectrum. For real pulses, only the positive part of the spectrum is shown

  • wmin (float) – Lowest frequency to show. Overrides zoom options. Must be given together with wmax.

  • wmax (float) – Highest frequency to show. Overrides zoom options. Must be given together with wmin.

  • spec_scale (float) – Factor by which to scale the amplitudes in the spectrum

  • spec_max (float) – Maximum amplitude in the spectrum, after spec_scale has been applied

  • freq_unit (str) – Unit in which to show the frequency axis in the spectrum. If not given, use the freq_unit attribute

  • mark_freqs (None, list(float), list((float, dict))) – Array of frequencies to mark in spectrum as vertical dashed lines. If list of tuples (float, dict), the float value is the frequency to mark, and the dict gives the keyword arguments that are passed to the matplotlib axvline method.

  • mark_freq_points (None, MarkerStyle) – Marker to be used to indicate the individual points in the spectrum.

The remaining figargs are passed to matplotlib.pyplot.figure to create a new figure if fig is None.

show(**kwargs)[source]

Show a plot of the pulse and its spectrum. All arguments will be passed to the plot method

show_pulse(**kwargs)[source]

Show a plot of the pulse amplitude; alias for show(show_spectrum=False). All other arguments will be passed to the show method

show_spectrum(zoom=True, freq_unit=None, **kwargs)[source]

Show a plot of the pulse spectrum; alias for show(show_pulse=False, zoom=zoom, freq_unit=freq_unit). All other arguments will be passed to the show method

qdyn.pulse.pulse_tgrid(T, nt, t0=0.0)[source]

Return a pulse time grid suitable for an equidistant time grid of the states between t0 and T with nt intervals. The values of the pulse are defined in the intervals of the time grid, so the pulse time grid will be shifted by dt/2 with respect to the time grid of the states. Also, the pulse time grid will have nt-1 points:

>>> print(", ".join([("%.2f" % t) for t in pulse_tgrid(1.5, nt=4)]))
0.25, 0.75, 1.25

The limits of the states time grid are defined as the starting and end points of the pulse, however:

>>> p = Pulse(tgrid=pulse_tgrid(1.5, 4), time_unit='ns', ampl_unit='MHz')
>>> p.t0
0
>>> p.T
1.5_ns
qdyn.pulse.tgrid_from_config(tgrid_dict, time_unit, pulse_grid=True)[source]

Extract the time grid from the given config file

>>> tgrid_dict = dict([('t_start', 0.0), ('t_stop', UnitFloat(10.0, 'ns')),
...                    ('dt', UnitFloat(20, 'ps')), ('fixed', True)])
>>> tgrid = tgrid_from_config(tgrid_dict, time_unit='ns')
>>> print("%.2f" % tgrid[0])
0.01
>>> print("%.2f" % tgrid[-1])
9.99
qdyn.pulse.carrier(t, time_unit, freq, freq_unit, weights=None, phases=None, complex=False)[source]

Create the “carrier” of the pulse as a weighted superposition of cosines at different frequencies.

Parameters
  • t (numpy.ndarray(float)) – Time value or time grid

  • time_unit (str) – Unit of t

  • freq (numpy.ndarray(float)) – Carrier frequency or frequencies

  • freq_unit (str) – Unit of freq

  • weights (numpy.ndarray) – If freq is an array, weights for the different frequencies. If not given, all weights are 1. The weights are normalized to sum to one. Any weight smaller than machine precision is assumed zero.

  • phases (numpy.ndarray) – If phases is an array, phase shift for each frequency component, in units of pi. If not given, all phases are 0.

  • complex (bool) – If True, oscillate in the complex plane

Returns

Depending on whether complex is True or False,

\[\begin{split}s(t) = \sum_j w_j * \cos(\omega_j * t + \phi_j) \\ s(t) = \sum_j w_j * \exp(i*(\omega_j * t + \phi_j))\end{split}\]

with \(\omega_j = 2 * \pi * f_j\), and frequency \(f_j\) where \(f_j\) is the j’th value in freq. The value of \(\phi_j\) is the j’th value in phases

signal is a scalar if t is a scalar, and and array if t is an array

Return type

numpy.ndarray(complex)

Notes

freq_unit can be Hz (GHz, MHz, etc), describing the frequency directly, or any energy unit, in which case the energy value E (given through the freq parameter) is converted to an actual frequency as

\[f = E / (\hbar * 2 * \pi)\]
qdyn.pulse.CRAB_carrier(t, time_unit, freq, freq_unit, a, b, normalize=False, complex=False)[source]

Construct a “carrier” based on the CRAB formula

\[E(t) = \sum_{n} (a_n \cos(\omega_n t) + b_n \cos(\omega_n t))\]

where \(a_n\) is the n’th element of a, \(b_n\) is the n’th element of b, and \(\omega_n\) is the n’th element of freq.

Parameters
  • t (numpy.ndarray) – time grid values

  • time_unit (str) – Unit of t

  • freq (numpy.ndarray) – Carrier frequency or frequencies

  • freq_unit (str) – Unit of freq

  • a (numpy.ndarray) – Coefficients for cosines

  • b (numpy.ndarray) – Coefficients for sines

  • normalize (bool) – If True, normalize the resulting carrier such that its values are in [-1,1]

  • complex (bool) –

    If True, oscillate in the complex plane

    \[E(t) = \sum_{n} (a_n - i b_n) \exp(i \omega_n t)\]

Notes

freq_unit can be Hz (GHz, MHz, etc), describing the frequency directly, or any energy unit, in which case the energy value E (given through the freq parameter) is converted to an actual frequency as

\[f = E / (\hbar * 2 * \pi)\]
qdyn.pulse.gaussian(t, t0, sigma)[source]

Return a Gaussian shape with peak amplitude 1.0

Parameters
Returns

Gaussian shape of same type as t

Return type

(float, numpy.ndarray)

qdyn.pulse.blackman(t, t_start, t_stop, a=0.16)[source]

Return a Blackman function between t_start and t_stop, see http://en.wikipedia.org/wiki/Window_function#Blackman_windows

A Blackman shape looks nearly identical to a Gaussian with a 6-sigma interval between start and stop Unlike the Gaussian, however, it will go exactly to zero at the edges. Thus, Blackman pulses are often preferable to Gaussians.

Parameters
  • t (float, numpy.ndarray) – Time point or time grid

  • t_start (float) – Starting point of Blackman shape

  • t_stop (float) – End point of Blackman shape

Returns

If t is a scalar, blackman_shape is the scalar value of the Blackman shape at t. If t is an array, blackman_shape is an array of same size as t, containing the values for the Blackman shape (zero before t_start and after t_stop)

Return type

(float, numpy.ndarray(float))

See also

numpy.blackman