quspin.tools.evolution.ExpmMultiplyParallel

class quspin.tools.evolution.ExpmMultiplyParallel(A, a=1.0, dtype=None, copy=False)[source]

Bases: object

Implements scipy.sparse.linalg.expm_multiply() for openmp.

Notes

  • this is a wrapper over custom c++ code.

  • the dtype input need not be the same dtype as A or a; however, it must be possible to cast the result of a*A to this dtype.

  • consider the special case of real-time evolution with a purely-imaginary Hamiltonian, in which case a=-1j*time and A are both complex-valued, while the resulting matrix exponential is real-valued: in such cases, one can use either one of

>>> expm_multiply_parallel( (1j*H.tocsr()).astype(np.float64), a=-1.0, dtype=np.float64)`

and

>>> expm_multiply_parallel( H.tocsr(), a=-1.0j, dtype=np.complex128)

The more efficient way to compute the matrix exponential in this case is to use a real-valued dtype.

Examples

This example shows how to construct the expm_multiply_parallel object.

Further code snippets can be found in the examples for the function methods of the class. The code snippet below initiates the class, and is required to run the example codes for the function methods.

 1from quspin.basis import spin_basis_1d  # bosonic Hilbert space
 2from quspin.tools.evolution import expm_multiply_parallel  # expm_multiply_parallel
 3import numpy as np  # general math functions
 4
 5#
 6L = 12  # syste size
 7# coupling strenghts
 8J = 1.0  # spin-spin coupling
 9h = 0.8945  # x-field strength
10g = 0.945  # z-field strength
11# create site-coupling lists
12J_zz = [[J, i, (i + 1) % L] for i in range(L)]  # PBC
13x_field = [[h, i] for i in range(L)]
14z_field = [[g, i] for i in range(L)]
15# create static and dynamic lists
16static = [["zz", J_zz], ["x", x_field], ["z", z_field]]
17dynamic = []
18# create spin-1/2 basis
19basis = spin_basis_1d(L, kblock=0, pblock=1)
20# set up Hamiltonian
21H = hamiltonian(static, dynamic, basis=basis, dtype=np.float64)
22# prealocate computation of matrix exponential
23expH = expm_multiply_parallel(H.tocsr(), a=0.2j)
24print(expH)
__init__(A, a=1.0, dtype=None, copy=False)[source]

Initializes expm_multiply_parallel.

Parameters:
A{array_like, scipy.sparse matrix}

The operator (matrix) whose exponential is to be calculated.

ascalar, optional

scalar value multiplying generator matrix \(A\) in matrix exponential: \(\mathrm{e}^{aA}\).

dtypenumpy.dtype, optional

data type specified for the total operator \(\mathrm{e}^{aA}\). Default is: numpy.result_type(A.dtype,min_scalar_type(a),float64).

copybool, optional

if True the matrix is copied otherwise the matrix is stored by reference.

Methods

__init__(A[, a, dtype, copy])

Initializes expm_multiply_parallel.

dot(v[, work_array, overwrite_v, tol])

Calculates the action of \(\mathrm{e}^{aA}\) on a vector \(v\).

set_a(a[, dtype])

Sets the value of the property a.

Attributes

A

csr_matrix to be exponentiated.

a

\(\mathrm{e}^{aA}\)

property A

csr_matrix to be exponentiated.

Type:

scipy.sparse.csr_matrix

property a

\(\mathrm{e}^{aA}\)

Type:

scalar

Type:

value multiplying generator matrix \(A\) in matrix exponential

dot(v, work_array=None, overwrite_v=False, tol=None)[source]

Calculates the action of \(\mathrm{e}^{aA}\) on a vector \(v\).

Parameters:
vcontiguous numpy.ndarray, 1d or 2d array

array to apply \(\mathrm{e}^{aA}\) on.

work_arraycontiguous numpy.ndarray, optional

array can be any shape but must contain 2*v.size contiguous elements. This array is used as temporary memory space for the underlying c-code. This saves extra memory allocation for function operations.

overwrite_vbool, optoinal

if set to True, the data in v is overwritten by the function. This saves extra memory allocation for the results.

tol: float, optoinal

tolerance value used to truncate Taylor expansion of matrix exponential.

Returns:
numpy.ndarray

result of \(\mathrm{e}^{aA}v\).

If overwrite_v = True the dunction returns v with the data overwritten, otherwise the result is stored in a new array.

Examples

 1##### compute expm_multiply applied on a state
 2_, psi = H.eigsh(k=1, which="SA")  # compute GS of H
 3psi = psi.squeeze().astype(
 4    np.complex128
 5)  # cast array type to complex double due to complex matrix exp
 6# construct c++ work array for speed
 7work_array = np.zeros((2 * len(psi),), dtype=psi.dtype)
 8print(H.expt_value(psi))  # measure energy of state |psi>
 9expH.dot(
10    psi, work_array=work_array, overwrite_v=True
11)  # compute action of matrix exponential on a state
12print(H.expt_value(psi))  # measure energy of state exp(aH)|psi>
set_a(a, dtype=None)[source]

Sets the value of the property a.

Parameters:
ascalar

new value of a.

dtypenumpy.dtype, optional

dtype specified for this operator. Default is: result_type(A.dtype,min_scalar_type(a),float64)

Examples

1##### change value of `a`
2print(expH.a)
3expH.set_a(0.3j)
4print(expH.a)