quimb.operator.builder¶
Given a single definition of hilbert space and set of terms, build out various representations of the operator.
This includes:
Sparse matrix in various formats (CSR, CSC, COO, etc.)
Dense matrix
Matrix product operator (MPO) for DMRG etc.
Dict of k-local terms (for use with PEPS etc.)
Coupling function, for use with VMC etc.
Attributes¶
Classes¶
Object for building operators with sparse structure. Specifically, |
Functions¶
|
|
|
|
|
Transform the terms in this operator by pre-prending pauli Z |
|
Simplify a sequence of operators acting on the same site. |
|
Simplify the given terms by combining operators acting on the same site, |
|
Decompose the given operator (specified as a label) into a sum of |
|
Decompose the given terms into a sum of Pauli strings. |
|
|
|
Create a sparse nested dictionary of how each term couples each |
|
Module Contents¶
- quimb.operator.builder._OPMAP¶
- quimb.operator.builder._OPCOMPLEX¶
- quimb.operator.builder.jordan_wigner_transform(terms, site_to_reg=None, reg_to_site=None)[source]¶
Transform the terms in this operator by pre-prending pauli Z strings to all creation and annihilation operators. This is always applied directly to the raw terms, so that the any fermionic ordering is respected. Any further transformations (e.g. simplification or pauli decomposition) should thus be applied after this transformation.
Note this doesn’t decompose +, - into (X + iY) / 2 and (X - iY) / 2, it just prepends Z strings. Call pauli_decompose after this to get the full decomposition.
The placement of the Z strings is defined by the ordering supplied by site_to_reg and reg_to_site, by default, it assumes the terms already are specified on a linear range of integers.
- Parameters:
terms (dict[term, coeff]) – The terms of the operator. Each term is a tuple of tuples, where each inner tuple is a pair of
(operator, site).site_to_reg (callable, optional) – A function that maps a site to a linear register index. If not provided, the sites are assumed to be linear integers already.
reg_to_site (callable, optional) – A function that maps a linear register index to a site. If not provided, the sites are assumed to be linear integers already.
- Returns:
terms_jordan_wigner – The transformed terms of the operator. Each term is a tuple of tuples, where each inner tuple is a pair of
(operator, site).- Return type:
dict[term, coeff]
- quimb.operator.builder.simplify_single_site_ops(coeff, ops)[source]¶
Simplify a sequence of operators acting on the same site.
- Parameters:
- Returns:
new_coeff (float or complex) – The simplified coefficient.
new_op (str) – The single, simplified operator that the sequence maps to, up to scalar multiplication.
Examples
>>> simplify_single_site_ops(1.0, ('+', 'z', 'z', 'z', 'z', '-')) (1.0, 'n')
>>> simplify_single_site_ops(1.0, ("x", "y", "z")) (-1j, 'I')
- quimb.operator.builder.simplify(terms, atol=1e-12, site_to_reg=None)[source]¶
Simplify the given terms by combining operators acting on the same site, putting them in canonical order, removing null-terms, and combining equivalent operator strings.
- Parameters:
terms (dict[term, coeff]) – The terms of the operator. Each term is a tuple of tuples, where each inner tuple is a pair of
(operator, site).atol (float, optional) – The absolute tolerance for considering coefficients after simplification to be null.
site_to_reg (callable, optional) – A function that maps a site to a linear register index. If not provided, the sites are assumed to be linear integers already. This is just used to sort the operators into canonical order.
- Returns:
terms_simplified – The simplified terms of the operator. Each term is a tuple of tuples, where each inner tuple is a pair of
(operator, site).- Return type:
dict[term, coeff]
- quimb.operator.builder.get_pauli_decomp(op, atol=1e-12, use_zx=False)[source]¶
Decompose the given operator (specified as a label) into a sum of Pauli components.
- Parameters:
- Returns:
The decomposition of the operator into Pauli components. Each tuple contains the coefficient and the Pauli operator label.
- Return type:
- quimb.operator.builder.pauli_decompose(terms, atol=1e-12, use_zx=False, site_to_reg=None)[source]¶
Decompose the given terms into a sum of Pauli strings.
- Parameters:
terms (dict[term, coeff]) – The terms of the operator. Each term is a tuple of tuples, where each inner tuple is a pair of
(operator, site).atol (float, optional) – The absolute tolerance for considering coefficients after decomposition to be null.
use_xz (bool, optional) – Whether to decompose in terms of the real ⴵ = ZX = iY instead of Y.
site_to_reg (callable, optional) – A function that maps a site to a linear register index. If not provided, the sites are assumed to be linear integers already. This is just used to sort the operators into canonical order.
- Returns:
terms_pauli_decomposed – The decomposed terms of the operator. Each term is a tuple of tuples, where each inner tuple is a pair of
(operator, site).- Return type:
dict[term, coeff]
- quimb.operator.builder.build_coupling_numba(term_store, site_to_reg, dtype=None)[source]¶
Create a sparse nested dictionary of how each term couples each local site configuration to which other local site configuration, and with what coefficient, suitable for use with numba.
- class quimb.operator.builder.SparseOperatorBuilder(terms=(), hilbert_space: quimb.operator.hilbertspace.HilbertSpace = None, dtype=None, jordan_wigner=False, pauli_decompose=False, atol=1e-12)[source]¶
Object for building operators with sparse structure. Specifically, a sum of terms, where each term is a product of operators, where each of these local operators acts on a single site and has at most one entry per row.
- Parameters:
terms (sequence, optional) – The terms to initialize the builder with.
add_termis simply called on each of these.hilbert_space (HilbertSpace) – The Hilbert space to build the operator in. If this is not supplied then a minimal Hilbert space will be constructed from the sites used, when required. The default symmetry and sector to build operators in is inherited from this Hilbert space, but can be overridden.
dtype (numpy.dtype or str, optional) – A default data type for the coefficients of the operator. If not provided, will be automatically determined at building time based on the terms in the operator. If the operator is complex, will be set to
np.complex128. If the operator is real, will be set tonp.float64. Individual building methods can override this.jordan_wigner (bool, optional) – Whether to apply the Jordan-Wigner transformation to the terms automatically when processing them. This prepends pauli Z strings to all creation and annihilation operators.
pauli_decompose (bool or "zx", optional) – Whether to apply the Pauli decomposition to the terms automatically when processing them. This decomposes all local operators into sums of Pauli operators. If “zx” is supplied, the decomposition is done in terms of the real ZX = iY operator instead of Y.
atol (float, optional) – The absolute tolerance for considering coefficients to be null when simplifying and decomposing terms.
- _sites_used¶
- _hilbert_space = None¶
- _terms_raw¶
- _terms_final = None¶
- _transform_jordan_wigner = False¶
- _transform_pauli_decompose = False¶
- _atol = 1e-12¶
- _dtype = None¶
- _coupling_maps¶
- _cache¶
- property sites_used¶
A tuple of the sorted coordinates/sites seen so far.
- property hilbert_space: quimb.operator.hilbertspace.HilbertSpace¶
The Hilbert space of the operator. Created from the sites seen so far if not supplied at construction.
- property nsites¶
The total number of coordinates/sites seen so far.
- property terms_raw¶
A tuple of the raw terms seen so far, as a mapping from operator strings to coefficients.
- _reset_caches()[source]¶
Reset any cached representations of the operator, used whenever the terms are modified in any way, and thus require reprocessing.
- _get_terms_final()[source]¶
Get the processed terms, applying any requested transformations, if not already done.
- property terms¶
A tuple of the, possibly transformed, terms seen so far.
- property nterms¶
The total number of terms seen so far.
- property locality¶
The locality of the operator, the maximum support of any term.
- property iscomplex¶
Whether the operator has complex terms.
- add_term(*coeff_ops)[source]¶
Add a term to the operator.
- Parameters:
coeff (float, optional) – The overall coefficient of the term.
ops (sequence of tuple[str, hashable]) –
The operators of the term, together with the sites they act on. Each term should be a pair of
(operator, site), whereoperatorcan be:'x','y','z': Pauli matrices'sx','sy','sz': spin operators (i.e. scaled Pauli matrices)'+','-': creation/annihilation operators'n','sn', or'h': number, symmetric number (n - 1/2) and hole (1 - n) operators.
And
siteis a hashable object that represents the site that the operator acts on. If this builder has an associated Hilbert space already, the site must be present in that Hilbert space, else a minimal Hilbert space will be constructed from the sites used, when required.
- jordan_wigner_transform(value=None)[source]¶
Toggle transforming the terms in this operator by pre-prending pauli Z strings to all creation and annihilation operators. This is always applied directly as the first processing step to the raw terms, so that the fermionic ordering is respected.
Note this doesn’t decompose +, - into (X + iY) / 2 and (X - iY) / 2, it just prepends Z strings. Use pauli_decompose to get the full decomposition.
The placement of the Z strings is defined by the ordering of the hilbert space, by default, the sorted order of the site labels.
- Parameters:
value (bool, optional) – Whether to apply the Jordan-Wigner transformation. If None (the default) then this method acts as toggle. Whereas supplying True or False explicitly sets or unsets the transformation.
- pauli_decompose(value=None, atol=None, use_zx=False)[source]¶
Transform the terms in this operator by decomposing them into Pauli strings.
- Parameters:
value (bool, optional) – Whether to apply the Pauli decomposition. If None (the default) then this method acts as toggle. Whereas supplying True or False explicitly sets or unsets the transformation.
- get_dtype(dtype=None)[source]¶
Calculate the numpy data type of the operator to use.
- Parameters:
dtype (numpy.dtype or str, optional) – The data type of the coefficients. If not provided, will be automatically determined based on the terms in the operator.
- Returns:
dtype – The data type of the coefficients.
- Return type:
- get_coupling_map(dtype=None)[source]¶
Build and cache the coupling map for the specified dtype.
- Parameters:
dtype (numpy.dtype, optional) – The data type of the coefficients. If not provided, will be automatically determined based on the terms in the operator.
- Returns:
coupling_map – The operator defined as tuple of flat arrays.
- Return type:
tuple[ndarray]
- flatconfig_coupling(flatconfig, dtype=None)[source]¶
Get an array of other configurations coupled to the given individual
flatconfigby this operator, and the corresponding coupling coefficients. This is for use with VMC for example.- Parameters:
flatconfig (array[np.uint8]) – The linear array of the configuration to get the coupling for.
dtype (numpy.dtype, optional) – The data type of the matrix. If not provided, will be automatically determined based on the terms in the operator.
- Returns:
coupled_flatconfigs (ndarray[np.uint8]) – Each distinct flatconfig coupled to
flatconfig.coeffs (ndarray[dtype]) – The corresponding coupling coefficients.
- config_coupling(config, dtype=None)[source]¶
Get a list of other configurations coupled to
configby this operator, and the corresponding coupling coefficients. This is for use with VMC for example.
- evaluate_exact_flatconfigs(fn_amplitude, progbar=False)[source]¶
Calculate the expectation value of this operator with respect to a wavefunction provided as a function with signature:
fn_amplitude(flatconfig: ndarray[np.uint8]) -> z: float | complex
- evaluate_exact_configs(fn_amplitude, progbar=False)[source]¶
Calculate the expectation value of this operator with respect to a wavefunction provided as a function with signature:
fn_amplitude(config: dict[site, int]) -> z: float | complex
- build_coo_data(sector=None, symmetry=None, dtype=None, parallel=False)[source]¶
Build the raw data for a sparse matrix in COO format, optionally in parallel.
- Parameters:
sector ({None, str, int, ((int, int), (int, int))}, optional) – The sector of the Hilbert space. If None, the default sector is used.
symmetry ({None, "Z2", "U1", "U1U1"}, optional) – The symmetry of the Hilbert space. If None, the default symmetry is used, or inferred from the supplied sector if possible.
dtype (numpy.dtype, optional) – The data type of the matrix. If not provided, will be automatically determined based on the terms in the operator.
parallel (bool or int, optional) – Whether to build the matrix in parallel (multi-threaded). If True, will use number of threads equal to the number of available CPU cores. If False, will use a single thread. If an integer is provided, it will be used as the number of threads to use.
- Returns:
data (array) – The data entries for the sparse matrix in COO format.
rows (array) – The row indices for the sparse matrix in COO format.
cols (array) – The column indices for the sparse matrix in COO format.
d (int) – The total number of basis states.
- build_sparse_matrix(sector=None, symmetry=None, dtype=None, stype='csr', parallel=False)[source]¶
Build a sparse matrix in the given format. Optionally in parallel.
- Parameters:
sector ({None, str, int, ((int, int), (int, int))}, optional) – The sector of the Hilbert space. If None, the default sector is used.
symmetry ({None, "Z2", "U1", "U1U1"}, optional) – The symmetry of the Hilbert space. If None, the default symmetry is used, or inferred from the supplied sector if possible.
dtype (numpy.dtype, optional) – The data type of the matrix. If not provided, will be automatically determined based on the terms in the operator.
stype (str, optional) – The sparse matrix format to use. Can be one of ‘coo’, ‘csr’, ‘csc’, ‘bsr’, ‘lil’, ‘dok’, or ‘dia’. Default is ‘csr’.
parallel (bool, optional) – Whether to build the matrix in parallel (multi-threaded).
- Return type:
scipy.sparse matrix
- build_dense(sector=None, symmetry=None, dtype=None, parallel=False)[source]¶
Get the dense (numpy.ndarray) matrix representation of this operator.
- Parameters:
sector ({None, str, int, ((int, int), (int, int))}, optional) – The sector of the Hilbert space. If None, the default sector is used.
symmetry ({None, "Z2", "U1", "U1U1"}, optional) – The symmetry of the Hilbert space. If None, the default symmetry is used, or inferred from the supplied sector if possible.
dtype (numpy.dtype, optional) – The data type of the matrix. If not provided, will be automatically determined based on the terms in the operator.
parallel (bool or int, optional) – Whether to build the matrix in parallel (multi-threaded). If True, will use number of threads equal to the number of available CPU cores. If False, will use a single thread. If an integer is provided, it will be used as the number of threads to use.
- Returns:
A – The dense matrix representation of this operator.
- Return type:
- matvec(x, out=None, sector=None, symmetry=None, dtype=None, parallel=False)[source]¶
Apply this operator lazily (i.e. without constructing a sparse matrix) to a vector. This uses less memory but is much slower.
- Parameters:
x (array) – The vector to apply the operator to.
out (array, optional) – An array to store the result in. If not provided, a new array will be created.
sector ({None, str, int, ((int, int), (int, int))}, optional) – The sector of the Hilbert space. If None, the default sector is used. The implicit size should match that of x.
symmetry ({None, "Z2", "U1", "U1U1"}, optional) – The symmetry of the Hilbert space. If None, the default symmetry is used, or inferred from the supplied sector if possible.
dtype (numpy.dtype, optional) – The data type of the matrix. If not provided, will be automatically set as the same as the input vector.
parallel (bool or int, optional) – Whether to apply the operator in parallel (multi-threaded). If True, will use number of threads equal to the number of available CPU cores. If False, will use a single thread. If an integer is provided, it will be used as the number of threads to use. Uses num_threads more memory but is faster.
- Returns:
out – The result of applying the operator to the vector.
- Return type:
- aslinearoperator(sector=None, symmetry=None, dtype=None, parallel=False)[source]¶
Get a scipy.sparse.linalg.LinearOperator for this operator. This is a lazy representation of the operator, which uses matvec to apply the operator to a vector. Less memory is required than constructing the full sparse matrix, but it is significantly slower.
Note currently the operator is assumed to be hermitian.
- Parameters:
sector ({None, str, int, ((int, int), (int, int))}, optional) – The sector of the Hilbert space. If None, the default sector is used.
symmetry ({None, "Z2", "U1", "U1U1"}, optional) – The symmetry of the Hilbert space. If None, the default symmetry is used, or inferred from the supplied sector if possible.
dtype (numpy.dtype, optional) – The data type of the matrix. If not provided, will be automatically determined based on the terms in the operator.
parallel (bool or int, optional) – Whether to apply the operator in parallel (multi-threaded). If True, will use number of threads equal to the number of available CPU cores. If False, will use a single thread. If an integer is provided, it will be used as the number of threads to use. Uses num_threads more memory but is faster.
- Returns:
Alo – The linear operator representation of this operator.
- Return type:
- build_local_terms(dtype=None)[source]¶
Get a dictionary of local terms, where each key is a sorted tuple of sites, and each value is the local matrix representation of the operator on those sites. For use with e.g. tensor network algorithms.
Note terms acting on the same sites are summed together and the size of each local matrix is exponential in the locality of that term.
- Returns:
Hk – The local terms.
- Return type:
dict[tuple[hashable], numpy.ndarray]
- build_local_ham(dtype=None)[source]¶
Get a LocalHamGen object for this operator.
- Parameters:
dtype (numpy.dtype, optional) – The data type of the matrix. If not provided, will be automatically determined based on the terms in the operator.
- Returns:
H – The local Hamiltonian representation of this operator.
- Return type:
- draw_state_machine(method='greedy', figsize='auto', G=None)[source]¶
Draw the fintie state machine for this operator, as if buildling the MPO.
- build_mpo(method='greedy', dtype=None, **mpo_opts)[source]¶
Build a matrix product operator (MPO) representation of this operator.
- Parameters:
method (str, optional) – The method to use for building the MPO. Currently only “greedy” is supported.
dtype (type, optional) – The data type of the MPO. If not supplied, will be chosen automatically based on the terms in the operator.
mpo_opts (keyword arguments) – Additional options to pass to the MPO constructor. See MatrixProductOperator for details.
- Returns:
mpo – The MPO representation of this operator.
- Return type: