Constructing time dependent systems

QobjEvo

class QobjEvo

A class for representing time-dependent quantum objects, such as quantum operators and states.

Importantly, QobjEvo instances are used to represent such time-dependent quantum objects when working with QuTiP solvers.

A QobjEvo instance may be constructed from one of the following:

  • a callable f(t: double, args: dict) -> Qobj that returns the value of the quantum object at time t.

  • a [Qobj, Coefficient] pair, where the Coefficient may be any item that coefficient can accept (e.g. a function, a numpy array of coefficient values, a string expression).

  • a Qobj (which creates a constant QobjEvo term).

  • a list of such callables, pairs or Qobjs.

  • a QobjEvo (in which case a copy is created, all other arguments are ignored except args which, if passed, replaces the existing arguments).

Parameters:
Q_objectcallable, list or Qobj

A specification of the time-depedent quantum object. See the paragraph above for a full description and the examples section below for examples.

argsdict, optional

A dictionary that contains the arguments for the coefficients. Arguments may be omitted if no function or string coefficients that require arguments are present.

tlistarray-like, optional

A list of times corresponding to the values of the coefficients supplied as numpy arrays. If no coefficients are supplied as numpy arrays, tlist may be omitted, otherwise it is required.

The times in tlist do not need to be equidistant, but must be sorted.

By default, a cubic spline interpolation will be used to interpolate the value of the (numpy array) coefficients at time t. If the coefficients are to be treated as step functions, pass the argument order=0 (see below).

orderint, default=3

Order of the spline interpolation that is to be used to interpolate the value of the (numpy array) coefficients at time t. 0 use previous or left value.

copybool, default=True

Whether to make a copy of the Qobj instances supplied in the Q_object parameter.

compressbool, default=True

Whether to compress the QobjEvo instance terms after the instance has been created.

This sums the constant terms in a single term and combines [Qobj, coefficient] pairs with the same Qobj into a single pair containing the sum of the coefficients.

See compress.

function_style{None, “pythonic”, “dict”, “auto”}

The style of function signature used by callables in Q_object. If style is None, the value of qutip.settings.core["function_coefficient_style"] is used. Otherwise the supplied value overrides the global setting.

boundary_conditions2-Tuple, str or None, optional

Boundary conditions for spline evaluation. Default value is None. Correspond to bc_type of scipy.interpolate.make_interp_spline. Refer to Scipy’s documentation for further details: https://docs.scipy.org/doc/scipy/reference/generated/scipy.interpolate.make_interp_spline.html

Attributes:
dimslist

List of dimensions keeping track of the tensor structure.

shape(int, int)

List of dimensions keeping track of the tensor structure.

typestr

Type of quantum object: ‘bra’, ‘ket’, ‘oper’, ‘operator-ket’, ‘operator-bra’, or ‘super’.

superrepstr

Representation used if type is ‘super’. One of ‘super’ (Liouville form) or ‘choi’ (Choi matrix with tr = dimension).

Examples

A QobjEvo constructed from a function:

def f(t, args):
    return qutip.qeye(N) * np.exp(args['w'] * t)

QobjEvo(f, args={'w': 1j})

For list based QobjEvo, the list must consist of Qobj or [Qobj, Coefficient] pairs:

QobjEvo([H0, [H1, coeff1], [H2, coeff2]], args=args)

The coefficients may be specified either using a Coefficient object or by a function, string, numpy array or any object that can be passed to the coefficient function. See the documentation of coefficient for a full description.

An example of a coefficient specified by a function:

def f1_t(t, args):
    return np.exp(-1j * t * args["w1"])

QobjEvo([[H1, f1_t]], args={"w1": 1.})

And of coefficients specified by string expressions:

H = QobjEvo(
    [H0, [H1, 'exp(-1j*w1*t)'], [H2, 'cos(w2*t)']],
    args={"w1": 1., "w2": 2.}
)

Coefficients maybe also be expressed as numpy arrays giving a list of the coefficient values:

tlist = np.logspace(-5, 0, 100)
H = QobjEvo(
    [H0, [H1, np.exp(-1j * tlist)], [H2, np.cos(2. * tlist)]],
    tlist=tlist
)

The coeffients array must have the same len as the tlist.

A QobjEvo may also be built using simple arithmetic operations combining Qobj with Coefficient, for example:

coeff = qutip.coefficient("exp(-1j*w1*t)", args={"w1": 1})
qevo = H0 + H1 * coeff
__call__()

Get the Qobj at t.

Parameters:
tfloat

Time at which the QobjEvo is to be evalued.

_argsdict [optional]

New arguments as a dict. Update args with arguments(new_args).

**kwargs

New arguments as a keywors. Update args with arguments(**new_args).

Notes

If both the positional _args and keywords are passed new values from both will be used. If a key is present with both, the _args dict value will take priority.

arguments(_args=None, **kwargs)

Update the arguments.

Parameters:
_argsdict [optional]

New arguments as a dict. Update args with arguments(new_args).

**kwargs

New arguments as a keywors. Update args with arguments(**new_args).

Notes

If both the positional _args and keywords are passed new values from both will be used. If a key is present with both, the _args dict value will take priority.

compress()

Look for redundance in the QobjEvo components:

Constant parts, (Qobj without Coefficient) will be summed. Pairs [Qobj, Coefficient] with the same Qobj are merged.

Example: [[sigmax(), f1], [sigmax(), f2]] -> [[sigmax(), f1+f2]]

The QobjEvo is transformed inplace.

Returns:
None
conj()

Get the element-wise conjugation of the quantum object.

copy()

Return a copy of this QobjEvo

dag()

Get the Hermitian adjoint of the quantum object.

dtype

Type of the data layers of the QobjEvo. When different data layers are used, we return the type of the sum of the parts.

expect(t, state, check_real=True)

Expectation value of this operator at time t with the state.

Parameters:
tfloat

Time of the operator to apply.

stateQobj

right matrix of the product

check_realbool (True)

Whether to convert the result to a real when the imaginary part is smaller than the real part by a dactor of settings.core['rtol'].

Returns:
expectfloat or complex

state.adjoint() @ self @ state if state is a ket. trace(self @ matrix) is state is an operator or operator-ket.

expect_data(t, state)

Expectation is defined as state.adjoint() @ self @ state if state is a vector, or state is an operator and self is a superoperator. If state is an operator and self is an operator, then expectation is trace(self @ matrix).

isbra

Indicates if the system represents a bra state.

isconstant

Does the system change depending on t

isket

Indicates if the system represents a ket state.

isoper

Indicates if the system represents an operator.

isoperbra

Indicates if the system represents a operator-bra state.

isoperket

Indicates if the system represents a operator-ket state.

issuper

Indicates if the system represents a superoperator.

linear_map(op_mapping, *, _skip_check=False)

Apply mapping to each Qobj contribution.

Example

QobjEvo([sigmax(), coeff]).linear_map(spre)

gives the same result has

QobjEvo([spre(sigmax()), coeff])

Parameters:
op_mapping: callable

Funtion to apply to each elements.

Returns:
QobjEvo

Modified object

Notes

Does not modify the coefficients, thus linear_map(conj) would not give the the conjugate of the QobjEvo. It’s only valid for linear transformations.

matmul(t, state)

Product of this operator at time t to the state. self(t) @ state

Parameters:
tfloat

Time of the operator to apply.

stateQobj

right matrix of the product

Returns:
productQobj

The result product as a Qobj

matmul_data(t, state, out=None)

Compute out += self(t) @ state

num_elements

Number of parts composing the system

tidyup(atol=1e-12)

Removes small elements from quantum object.

to(data_type)

Convert the underlying data store of all component into the desired storage representation.

The different storage representations available are the “data-layer types”. By default, these are Dense, Dia and CSR, which respectively construct a dense matrix, diagonal sparse matrixand a compressed sparse row one.

The QobjEvo is transformed inplace.

Parameters:
data_typetype

The data-layer type that the data of this Qobj should be converted to.

Returns:
None
to_list()

Restore the QobjEvo to a list form.

Returns:
list_qevo: list

The QobjEvo as a list, element are either Qobj for constant parts, [Qobj, Coefficient] for coefficient based term. The original format of the Coefficient is not restored. Lastly if the original QobjEvo is constructed with a function returning a Qobj, the term is returned as a pair of the original function and args (dict).

trans()

Transpose of the quantum object

Coefficient

coefficient(
base: CoefficientLike,
*,
tlist: ArrayLike = None,
args: dict = {},
args_ctypes: dict = {},
order: int = 3,
compile_opt: dict = None,
function_style: str = None,
boundary_conditions: tuple | str = None,
**kwargs,
)[source]

Build Coefficient for time dependent systems:

` QobjEvo = Qobj + Qobj * Coefficient + Qobj * Coefficient + ... `

The coefficients can be a function, a string or a numpy array. Other packages may add support for other kind of coefficients.

For function based coefficients, the function signature must be either:

  • f(t, ...) where the other arguments are supplied as ordinary “pythonic” arguments (e.g. f(t, w, a=5))

  • f(t, args) where the arguments are supplied in a “dict” named args

By default the signature style is controlled by the qutip.settings.core["function_coefficient_style"] setting, but it may be overriden here by specifying either function_style="pythonic" or function_style="dict".

Examples:

  • pythonic style function signature:

    def f1_t(t, w):
        return np.exp(-1j * t * w)
    
    coeff1 = coefficient(f1_t, args={"w": 1.})
    
  • dict style function signature:

    def f2_t(t, args):
        return np.exp(-1j * t * args["w"])
    
    coeff2 = coefficient(f2_t, args={"w": 1.})
    

For string based coeffients, the string must be a compilable python code resulting in a complex. The following symbols are defined:

sin, cos, tan, asin, acos, atan, pi, sinh, cosh, tanh, asinh, acosh, atanh, exp, log, log10, erf, zerf, sqrt, real, imag, conj, abs, norm, arg, proj, numpy as np, scipy.special as spe (python interface) and cython_special (scipy cython interface)

Examples:

coeff = coefficient('exp(-1j*w1*t)', args={"w1":1.})

‘args’ is needed for string coefficient at compilation. It is a dict of (name:object). The keys must be a valid variables string.

Compilation options can be passed as “compile_opt=CompilationOptions(…)”.

For numpy array format, the array must be an 1d of dtype float or complex. A list of times (float64) at which the coeffients must be given (tlist). The coeffients array must have the same len as the tlist. The time of the tlist do not need to be equidistant, but must be sorted. By default, a cubic spline interpolation will be used to compute the coefficient at time t. The keyword order sets the order of the interpolation. When order = 0, the interpolation is step function that evaluates to the most recent value.

Examples:

tlist = np.logspace(-5,0,100)
H = QobjEvo(np.exp(-1j*tlist), tlist=tlist)

scipy.interpolate’s CubicSpline, PPoly and Bspline are also converted to interpolated coefficients (the same kind of coefficient created from ndarray). Other interpolation methods from scipy are converted to a function-based coefficient (the same kind of coefficient created from callables).

Parameters:
baseobject

Base object to make into a Coefficient.

argsdict, optional

Dictionary of arguments to pass to the function or string coefficient.

orderint, default=3

Order of the spline for array based coefficient.

tlistiterable, optional

Times for each element of an array based coefficient.

function_stylestr {“dict”, “pythonic”, None}, optional

Function signature of function based coefficients.

args_ctypesdict, optional

C type for the args when compiling array based coefficients.

compile_optCompilationOptions, optional

Sets of options for the compilation of string based coefficients.

boundary_conditions: 2-tupule, str or None, optional

Specify boundary conditions for spline interpolation.

**kwargs

Extra arguments to pass the the coefficients.

CompilationOptions

class CompilationOptions(**options)[source]

Options that control compilation of string based coefficient to Cython.

These options can be set globaly:

settings.compile["compiler_flags"] = "-O1"

In a with block:

with CompilationOptions(use_cython=False):

Or as an instance:

coefficient(coeff, compile_opt=CompilationOptions(recompile=True))

Compilation options:

use_cython: bool

Whether to compile strings as cython code or use python’s exec.

recompilebool

Do not use previously made files but build a new one.

try_parse: bool [True]

Whether to try parsing the string for reuse and static typing.

static_typesbool [True]

Whether to use C types for constant and args.

accept_intNone, bool

Whether to use the type int for integer constants and args or upgrade it to float or complex. If None, it will only use int when subscription is found in the code.

accept_floatbool

Whether to use the type float or upgrade them to complex.

compiler_flagsstr

Flags to pass to the compiler, ex: “-Wall -O3”… Flags not matching your comiler and OS may cause compilation to fail. Use “recompile=True”, when trying to if the string pattern was previously used.

link_flagsstr

Libraries to link to pass to the compiler. They can not be used to add function to the string coefficient.

extra_importstr

Cython code to add at the head of the file. Can be used to add extra import or cimport code, ex: extra_import=”from scipy.linalg import det” extra_import=”from qutip.core.data cimport CSR”

clean_on_errorbool [True]

When writing a cython file that cannot be imported, erase it.

build_dir: str [None]

cythonize’s build_dir.