CCE Calculators
Documentation for the calculator objects called by Simulator
object.
Base class
- class RunObject(timespace, clusters, bath, magnetic_field, center=None, pulses=None, nbstates=None, seed=None, masked=True, parallel=False, direct=False, parallel_states=False, **kwargs)
Abstract class of the CCE simulation runner.
Implements cluster correlation expansion, interlaced averaging, and sampling over random bath states. Requires definition of the following methods, from which the kernel will be automatically created:
.generate_hamiltonian(self)
method which, using the attributes of theself
object, computes cluster hamiltonian stored inself.cluster_hamiltonian
..compute_result(self)
method which, using the attributes of theself
, computes the resulting quantity for the given cluster.
Alternatively, user can define the kernel manually. Then the following methods have to be overridden:
.kernel(self, cluster, *args, **kwargs)
method which takes indexes of the bath spins in the given cluster as a first positional argument. This method is required for usual CCE runs..interlaced_kernel(self, cluster, supercluster, *args, **kwargs)
method which takes indexes of the bath spins in the given cluster as a first positional argument, indexes of the supercluster as a second positional argument. This method is required for interlaced CCE runs.
- Parameters
timespace (ndarray with shape (t, )) – Time delay values at which to compute propagators.
clusters (dict) – Clusters included in different CCE orders of structure
{int order: ndarray([[i,j],[i,j]])}
.bath (BathArray with shape (n,)) – Array of n bath spins.
magnetic_field (ndarray) – Magnetic field of type
magnetic_field = np.array([Bx, By, Bz])
.alpha (int or ndarray with shape (2s+1, )) – \(\ket{0}\) state of the qubit in \(S_z\) basis or the index of eigenstate to be used as one.
beta (int or ndarray with shape (2s+1, )) – \(\ket{1}\) state of the qubit in \(S_z\) basis or the index of the eigenstate to be used as one.
state (ndarray with shape (2s+1, )) – Initial state of the central spin, used in gCCE and noise autocorrelation calculations. Defaults to \(\frac{1}{N}(\ket{0} + \ket{1})\) if not set OR if alpha and beta are provided as indexes.
spin (float) – Value of the central spin.
zfs (ndarray with shape (3,3)) – Zero Field Splitting tensor of the central spin.
gyro (float or ndarray with shape (3, 3)) –
Gyromagnetic ratio of the central spin
OR
tensor corresponding to interaction between magnetic field and central spin.
nbstates (int) – Number of random bath states to sample over in bath state sampling runs.
bath_state (ndarray) – Array of bath states in any accepted format.
seed (int) – Seed for the random number generator in bath states sampling.
masked (bool) – True if mask numerically unstable points (with result > result[0]) in the sampling over bath states False if not. Default True.
projected_bath_state (ndarray with shape (n,)) – Array with z-projections of the bath spins states. Overridden in runs with random bath state sampling.
parallel (bool) – True if parallelize calculation of cluster contributions over different mpi processes. Default False.
direct (bool) – True if use direct approach in run (requires way more memory but might be more numerically stable). False if use memory efficient approach. Default False.
parallel_states (bool) – True if use MPI to parallelize the calculations of density matrix for each random bath state.
**kwargs – Additional keyword arguments to be set as the attributes of the given object.
- result_operator(b, /)
Operator which will combine the result of expansion,.
Default:
operator.imul
.
- contribution_operator(b, /)
Operator which will combine multiple contributions of the same cluster in the optimized approach.
Default:
operator.ipow
.
- removal_operator(b, /)
Operator which will remove subcluster contribution from the given cluster contribution. First argument cluster contribution, second - subcluster contribution.
Defalut:
operator.itruediv
.
- addition_operator(axis=None, dtype=None, out=None, keepdims=<no value>, initial=<no value>, where=<no value>)
Group operation which will combine contributions from the different clusters into one contribution in the direct approach.
Default:
numpy.prod
.
- nbstates
Number of random bath states to sample over in bath state sampling runs.
- Type
int
- timespace
Time points at which result will be computed.
- Type
ndarray with shape (t, )
- clusters
Clusters included in different CCE orders of structure
{int order: ndarray([[i,j],[i,j]])}
.- Type
dict
- bath
Array of n bath spins.
- Type
BathArray with shape (n,)
- center
Properties of the central spin.
- Type
- magnetic_field
Magnetic field of type
magnetic_field = np.array([Bx, By, Bz])
, or a function that takes position as an argument.- Type
ndarray or callable
- parallel
True if parallelize calculation of cluster contributions over different mpi processes. Default False.
- Type
bool
- parallel_states
True if use MPI to parallelize the calculations of density matrix for each random bath state.
- Type
bool
- direct
True if use direct approach in run (requires way more memory but might be more numerically stable). False if use memory efficient approach. Default False.
- Type
bool
- seed
Seed for the random number generator in bath states sampling.
- Type
int
- masked
True if mask numerically unstable points (with result > result[0]) in the sampling over bath states False if not. Default True.
- Type
bool
- hamiltonian
Full cluster Hamiltonian.
- Type
ndarray
- has_states
Whether there are states provided in the bath during the run.
- Type
bool
- initial_states_mask
Bool array of the states, initially present in the bath.
- Type
ndarray
- projected_states
Array of \(S_z\) projections of the bath spins after each control pulse, involving bath spins.
- Type
ndarray
- base_hamiltonian
Hamiltonian of the given cluster without mean field additions. In conventional CCE, also excludes additions from central spins.
- Type
- result
Result of the calculation.
- Type
ndarray
- delays
List with delays before each pulse or None if equispaced. Generated by
.generate_pulses
method.- Type
list or None
- rotations
List with matrix representations of the rotation from each pulse. Generated by
.generate_pulses
method.- Type
list
- preprocess()
Method which will be called before cluster-expanded run.
- postprocess()
Method which will be called after cluster-expanded run.
- kernel(cluster, *args, **kwargs)
Central kernel that will be called in the cluster-expanded calculations.
- Parameters
cluster (ndarray) – Indexes of the bath spins in the given cluster.
*args – Positional arguments of the kernel.
**kwargs – Keyword arguments of the kernel.
- Returns
Results of the calculations.
- Return type
ndarray
- run(*args, **kwargs)
Method that runs cluster-expanded single calculation.
- Parameters
*args – Positional arguments of the kernel.
**kwargs – Keyword arguments of the kernel.
- Returns
Results of the calculations.
- Return type
ndarray
- sampling_run(*args, **kwargs)
Method that runs bath sampling calculations.
- Parameters
*args – Positional arguments of the kernel.
**kwargs – Keyword arguments of the kernel.
- Returns
Results of the calculations.
- Return type
ndarray
- interlaced_kernel(cluster, supercluster, *args, **kwargs)
Central kernel that will be called in the cluster-expanded calculations with interlaced averaging of bath spin states.
- Parameters
cluster (ndarray) – Indexes of the bath spins in the given cluster.
supercluster (ndarray) – Indexes of the bath spins in the supercluster of the given cluster. Supercluster is the union of all clusters in
.clusters
attribute, for which given cluster is a subset.*args – Positional arguments of the kernel.
**kwargs – Keyword arguments of the kernel.
- Returns
Results of the calculations.
- Return type
ndarray
- interlaced_run(*args, **kwargs)
Method that runs cluster-expanded single calculation with interlaced averaging of bath spin states.
- Parameters
*args – Positional arguments of the interlaced kernel.
**kwargs – Keyword arguments of the interlaced kernel.
- Returns
Results of the calculations.
- Return type
ndarray
- sampling_interlaced_run(*args, **kwargs)
Method that runs bath sampling calculations with interlaced averaging of bath spin states.
- Parameters
*args – Positional arguments of the interlaced kernel.
**kwargs – Keyword arguments of the interlaced kernel.
- Returns
Results of the calculations.
- Return type
ndarray
- classmethod from_simulator(sim, **kwargs)
Class method to generate
RunObject
from the properties ofSimulator
object.
- generate_supercluser_states(supercluster)
Helper function to generate all possible pure states of the given supercluster.
- Parameters
supercluster (ndarray with shape (n, )) – Indexes of the bath spins in the supercluster.
- Yields
ndarray with shape (n, ) – Pure state of the given supercluster.
- generate_pulses()
Generate list of matrix representations of the rotations, induced by the sequence of the pulses.
- Returns
tuple containing:
list or None: List with delays before each pulse or None if equispaced.
list: List with matrix representations of the rotation from each pulse.
- Return type
tuple
- get_hamiltonian_variable_bath_state(index=0)
Generate Hamiltonian in case of the complicated pulse sequence.
- Parameters
index (int) – Index of the flips of spin states.
- Returns
Hamiltonian with mean field additions from the given set of projected states.
- Return type
ndarray with shape (n, n)
- from_sigma(sigma, i, dims)
Generate spin vector from dictionary with spin matrices.
- Parameters
sigma (dict) – Dictionary, which contains spin matrices of form
{'x': Sx, 'y': Sy, 'z': Sz}
.i (int) – Index of the spin in the order of
dims
.dims (ndarray with shape (N,)) – Dimensions of the spins in the given cluster.
- Returns
Spin vector in a full Hilbert space.
- Return type
ndarray with shape (3, n, n)
- generate_rotated_projected_states(bath, pulses)
Generate projected states after each control pulse, involving bath spins.
- Parameters
bath (BathArray with shape (n, )) – Array of bath spins.
pulses (Sequence) – Sequence of pulses.
- Returns
Array of \(S_z\) projections of bath spin states after each pulse, involving bath spins. Each \(i\)-th column is projections before the \(i\)-th pulse involving bath spins.
- Return type
ndarray with shape (n, x)
- pulse_bath_rotation(pulse, vectors)
Generate rotation of the bath spins from the given pulse.
- Parameters
pulse (Pulse) – Control pulse.
vectors (ndarray with shape (n, 3, N, N) – Array of spin vectors.
- Returns
Matrix representation of the spin rotation.
- Return type
ndarray with shape (x, x)
- simple_propagator(timespace, hamiltonian)
Generate a simple propagator \(U=\exp[-\frac{i}{\hbar} \hat H]\) from the Hamiltonian.
- Parameters
timespace (ndarray with shape (n, )) – Time points at which to evaluate the propagator.
hamiltonian (ndarray with shape (N, N)) – Hamiltonian of the system.
- Returns
Propagators, evaluated at each timepoint.
- Return type
ndarray with shape (n, N, N)
- from_central_state(dimensions, central_state)
Generate density matrix of the system if all spins apart from central spin are in completely mixed state.
- Parameters
dimensions (ndarray with shape (n,)) – Array of the dimensions of the spins in the cluster.
central_state (ndarray with shape (x,)) – Density matrix of central spins.
- Returns
Density matrix for the whole cluster.
- Return type
ndarray with shape (N, N)
- from_none(dimensions)
Generate density matrix of the systems if all spins are in completely mixed state. :param dimensions: Array of the dimensions of the spins in the cluster. :type dimensions: ndarray with shape (n,)
- Returns
Density matrix for the whole cluster.
- Return type
ndarray with shape (N, N)
- from_states(states)
Generate density matrix of the systems if all spins are in pure states. :param states: Array of the pure spin states. :type states: array-like
- Returns
Spin vector for the whole cluster.
- Return type
ndarray with shape (N, N)
- combine_cluster_central(cluster_state, central_state)
Combine bath spin states and the state of central spin. :param cluster_state: State vector or density matrix of the bath spins. :type cluster_state: ndarray with shape (n,) or (n, n) :param central_state: State vector or density matrix of the central spins. :type central_state: ndarray with shape (m,) or (m, m)
- Returns
State vector or density matrix of the full system.
- Return type
ndarray with shape (mn, ) or (mn, mn)
- rand_state(d)
Generate random state of the spin.
- Parameters
d (int) – Dimensions of the spin.
- Returns
Density matrix of the random state.
- Return type
ndarray with shape (d, d)
- generate_initial_state(dimensions, states=None, central_state=None)
Generate initial state of the cluster.
- Parameters
dimensions (ndarray with shape (n, )) – Dimensions of all spins in the cluster.
states (BathState, optional) – States of the bath spins. If None, assumes completely random state.
central_state (ndarray) – State of the central spin. If None, assumes that no central spin is present in the Hilbert space of the cluster.
- Returns
State vector or density matrix of the cluster.
- Return type
ndarray with shape (N,) or (N, N)
Conventional CCE
- simple_propagators(delays, hamiltonian_alpha, hamiltonian_beta)
Generate two simple propagators \(U=\exp[-\frac{i}{\hbar} \hat H]\) from the Hamiltonians, conditioned on two qubit levels.
- Parameters
delays (ndarray with shape (n, )) – Time points at which to evaluate the propagator.
hamiltonian_alpha (ndarray with shape (N, N)) – Hamiltonian of the bath spins with qubit in alpha state.
hamiltonian_beta (ndarray with shape (N, N)) – Hamiltonian of the bath spins with qubit in beta state.
- Returns
ndarray with shape (n, N, N): Matrix representation of the propagator conditioned on the alpha qubit state for each time point.
ndarray with shape (n, N, N): Matrix representation of the propagator conditioned on the beta qubit state for each time point.
- Return type
tuple
- propagate_propagators(v0, v1, number)
From two simple propagators and number of pulses in CPMG sequence generate two full propagators. :param v0: Propagator conditioned on the alpha qubit state for each time point. :type v0: ndarray with shape (n, N, N) :param v1: Propagator conditioned on the beta qubit state for each time point. :type v1: ndarray with shape (n, N, N) :param number: Number of pulses. :type number: int
- Returns
ndarray with shape (n, N, N): Matrix representation of the propagator conditioned on the alpha qubit state for each time point.
ndarray with shape (n, N, N): Matrix representation of the propagator conditioned on the beta qubit state for each time point.
- Return type
tuple
- class CCE(*args, as_delay=False, second_order=False, level_confidence=0.95, **kwargs)
Class for running conventional CCE simulations.
Note
Subclass of the
RunObject
abstract class.- Parameters
*args – Positional arguments of the
RunObject
.pulses (int or Sequence) – number of pulses in CPMG sequence or instance of Sequence object. For now, only CPMG sequences are supported in conventional CCE simulations.
as_delay (bool) – True if time points are delay between pulses, False if time points are total time.
second_order (bool) – True if add second order perturbation theory correction to the cluster Hamiltonian. If set to True sets the qubit states as eigenstates of central spin Hamiltonian from the following procedure. If qubit states are provided as vectors in \(S_z\) basis, for each qubit state compute the fidelity of the qubit state and all eigenstates of the central spin and chose the one with fidelity higher than
level_confidence
. If such state is not found, raises an error.level_confidence (float) – Maximum fidelity of the qubit state to be considered eigenstate of the central spin Hamiltonian. Default 0.95.
**kwargs – Keyword arguments of the
RunObject
.
- pulses
If input Sequence contains only pi pulses at even delay, stores number of pulses. Otherwise stores full
Sequence
.- Type
int or Sequence
- as_delay
True if time points are delay between pulses, False if time points are total time.
- Type
bool
- second_order
True if add second order perturbation theory correction to the cluster hamiltonian.
- Type
bool
- level_confidence
Maximum fidelity of the qubit state to be considered eigenstate of the central spin hamiltonian.
- Type
float
- energy_alpha
Eigen energy of the alpha state in the central spin Hamiltonian.
- Type
float
- energy_beta
Eigen energy of the beta state in the central spin Hamiltonian.
- Type
float
- energies
All eigen energies of the central spin Hamiltonian.
- Type
ndarray with shape (2s+1, )
- projections_alpha_all
Array of vectors with spin operator matrix elements of type \([\bra{0}\hat S_x\ket{i}, \bra{0}\hat S_y\ket{i}, \bra{0}\hat S_z\ket{i}]\), where \(\ket{0}\) is the alpha qubit state, \(\ket{i}\) are all eigenstates of the central spin hamiltonian.
- Type
ndarray with shape (2s+1, 3)
- projections_beta_all
Array of vectors with spin operator matrix elements of type \([\bra{1}\hat S_x\ket{i}, \bra{1}\hat S_y\ket{i}, \bra{1}\hat S_z\ket{i}]\), where \(\ket{1}\) is the beta qubit state, \(\ket{i}\) are all eigenstates of the central spin hamiltonian.
- Type
ndarray with shape (2s+1, 3)
- projections_alpha
Vector with spin operator matrix elements of type \([\bra{0}\hat S_x\ket{0}, \bra{0}\hat S_y\ket{0}, \bra{0}\hat S_z\ket{0}]\), where \(\ket{0}\) is the alpha qubit state
- Type
ndarray with shape (3,)
- projections_beta
Vectors with spin operator matrix elements of type \([\bra{1}\hat S_x\ket{1}, \bra{1}\hat S_y\ket{1}, \bra{1}\hat S_z\ket{1}]\), where \(\ket{1}\) is the beta qubit state.
- Type
ndarray with shape (3,)
- use_pulses
True if use full
Sequence
. False if use only number of pulses.- Type
bool
- preprocess()
Method which will be called before cluster-expanded run.
- postprocess()
Method which will be called after cluster-expanded run.
- generate_hamiltonian()
Using the attributes of the
self
object, compute the two projected cluster hamiltonians.- Returns
Tuple containing:
Hamiltonian: Cluster hamiltonian when qubit in the alpha state.
Hamiltonian: Cluster hamiltonian when qubit in the alpha state.
- Return type
tuple
- compute_result()
Using the attributes of the
self
object, compute the coherence function as overlap in the bath evolution.- Returns
Computed coherence.
- Return type
ndarray
- propagators()
Generate two propagators, conditioned on the qubit state.
- Returns
tuple containing:
ndarray with shape (t, n, n): Matrix representation of the propagator conditioned on the alpha qubit state for each time point.
ndarray with shape (t, n, n): Matrix representation of the propagator conditioned on the beta qubit state for each time point.
- Return type
tuple
Generalized CCE
- rotation_propagator(u, rotations)
Generate the propagator from the simple propagator and set of \(2 au\) equispaced rotation operators.
Note
While the spacing between rotation operators is assumed to be \(2 au\), the spacing before and after the first and the last rotation respectively is assumed to be :math:` au`.
- Parameters
u (ndarray with shape (n, N, N)) – Simple propagator.
rotations (ndarray with shape (x, N, N)) – Array of rotation operators.
- Returns
Full propagator.
- Return type
ndarray with shape (n, N, N)
- class gCCE(*args, i=None, j=None, as_delay=False, fulldm=False, normalized=True, **kwargs)
Class for running generalized CCE simulations.
Note
Subclass of the
RunObject
abstract class.- Parameters
*args – Positional arguments of the
RunObject
.pulses (Sequence) – Sequence object, containing series of pulses, applied to the system.
as_delay (bool) – True if time points are delay between pulses, False if time points are total time.
fulldm (bool) – True if return full density matrix. Default False.
**kwargs – Keyword arguments of the
RunObject
.
- as_delay
True if time points are delay between pulses, False if time points are total time.
- Type
bool
- dm0
Initial density matrix of the central spin.
- Type
ndarray with shape (2s+1, 2s+1)
- normalization
Coherence at time 0.
- Type
float
- zero_cluster
Coherence computed for the isolated central spin.
- Type
ndarray with shape (n,)
- fulldm
True if return full density matrix.
- Type
bool
- preprocess()
Method which will be called before cluster-expanded run.
- process_dm(density_matrix)
Obtain the result from the density matrices.
- Parameters
density_matrix (ndarray with shape (n, N, N)) – Array of the density matrices.
- Returns
Depending on the parameters, returns the off diagonal element of the density matrix or full matrix.
- Return type
ndarray
- postprocess()
Method which will be called after cluster-expanded run.
- generate_hamiltonian()
Using the attributes of the
self
object, compute the cluster hamiltonian including the central spin.- Returns
Cluster hamiltonian.
- Return type
- compute_result()
Using the attributes of the
self
object, compute the coherence function of the central spin.- Returns
Computed coherence.
- Return type
ndarray
- propagator()
Function to compute time propagator U.
- Returns
Array of propagators, evaluated at each time point in
self.timespace
.- Return type
ndarray with shape (t, n, n)
Noise Autocorrelation
- correlation_it_j0(operator_i, operator_j, dm0_expanded, U)
Function to compute correlation function of the operator i at time t and operator j at time 0
- Parameters
operator_i (ndarray with shape (n, n)) – Matrix representation of operator i.
operator_j (ndarray with shape (n, n)) – Matrix representation of operator j.
dm0_expanded (ndarray with shape (n, n)) – Initial density matrix of the cluster.
U (ndarray with shape (t, n, n)) – Time evolution propagator, evaluated over t time points.
- Returns
Autocorrelation of the z-noise at each time point.
- Return type
ndarray with shape (t,)
- compute_correlations(nspin, dm0_expanded, U, central_spin=None)
Function to compute correlations for the given cluster, given time propagator U.
- Parameters
nspin (BathArray) – BathArray of the given cluster of bath spins.
dm0_expanded (ndarray with shape (n, n)) – Initial density matrix of the cluster.
U (ndarray with shape (t, n, n)) – Time evolution propagator, evaluated over t time points.
central_spin (CenterArray) – Array of central spins.
- Returns
correlation of the Overhauser field, induced by the given cluster at each time point.
- Return type
ndarray with shape (t,)
- class gCCENoise(*args, **kwargs)
Class for running generalized CCE simulations of the noise autocorrelation function.
Note
Subclass of the
RunObject
abstract class.- Parameters
*args – Positional arguments of the
RunObject
.**kwargs – Keyword arguments of the
RunObject
.
- result_operator(b, /)
Overridden operator which will combine the result of expansion:
operator.iadd
.
- contribution_operator(b, /)
Overridden operator which will combine multiple contributions of the same cluster in the optimized approach:
operator.imul
.
- removal_operator(b, /)
Overridden operator which remove subcluster contribution from the given cluster contribution:
operator.isub
.
- addition_operator(axis=None, dtype=None, out=None, keepdims=<no value>, initial=<no value>, where=<no value>)
Overridden group operation which will combine contributions from the different clusters into one contribution in the direct approach:
numpy.sum
.
- preprocess()
Method which will be called before cluster-expanded run.
- postprocess()
Method which will be called after cluster-expanded run.
- generate_hamiltonian()
Using the attributes of the
self
object, compute the cluster hamiltonian including the central spin.- Returns
Cluster hamiltonian.
- Return type
- compute_result()
Using the attributes of the
self
object, compute autocorrelation function of the noise from bath spins in the given cluster.- Returns
Computed autocorrelation function.
- Return type
ndarray
- class CCENoise(*args, **kwargs)
Class for running conventional CCE simulations of the noise autocorrelation function.
Note
Subclass of the
RunObject
abstract class.Warning
In general, for calculations of the autocorrelation function, better results are achieved with generalized CCE, which accounts for the evolution of the entangled state of the central spin.
Second order couplings between nuclear spins are not implemented.
- Parameters
*args – Positional arguments of the
RunObject
.**kwargs – Keyword arguments of the
RunObject
.
- result_operator(b, /)
Overridden operator which will combine the result of expansion:
operator.iadd
.
- contribution_operator(b, /)
Overridden operator which will combine multiple contributions of the same cluster in the optimized approach:
operator.imul
.
- removal_operator(b, /)
Overridden operator which remove subcluster contribution from the given cluster contribution:
operator.isub
.
- addition_operator(axis=None, dtype=None, out=None, keepdims=<no value>, initial=<no value>, where=<no value>)
Overridden group operation which will combine contributions from the different clusters into one contribution in the direct approach:
numpy.sum
.
- preprocess()
Method which will be called before cluster-expanded run.
- postprocess()
Method which will be called after cluster-expanded run.
- generate_hamiltonian()
Using the attributes of the
self
object, compute the projected cluster hamiltonian, averaged for two qubit states.- Returns
Cluster hamiltonian.
- Return type
- compute_result()
Using the attributes of the
self
object, compute autocorrelation function of the noise from bath spins in the given cluster.- Returns
Computed autocorrelation function.
- Return type
ndarray
Cluster-correlation Expansion Decorators
The way we find cluster in the code.
- generate_clusters(bath, r_dipole, order, r_inner=0, ignore=None, strong=False, nclusters=None)
Generate clusters for the bath spins.
- Parameters
bath (BathArray) – Array of bath spins.
r_dipole (float) – Maximum connectivity distance.
order (int) – Maximum size of the clusters to find.
r_inner (float) – Minimum connectivity distance.
ignore (list or str, optional) – If not None, includes the names of bath spins which are ignored in the cluster generation.
strong (bool) – Whether to find only completely interconnected clusters (default False).
nclusters (dict) –
Dictionary which contain maximum number of clusters of the given size. Has the form
n_clusters = {order: number}
, whereorder
is the size of the cluster,number
is the maximum number of clusters with this size.If provided, sorts the clusters by the strength of cluster interaction, equal to the lowest pairwise interaction in the cluster. Then the strongest
number
of clusters is taken.
- Returns
Dictionary with keys corresponding to size of the cluster, and value corresponds to ndarray of shape (matrix, N). Here matrix is the number of clusters of given size, N is the size of the cluster. Each row contains indexes of the bath spins included in the given cluster.
- Return type
dict
- make_graph(bath, r_dipole, r_inner=0, ignore=None, max_size=5000)
Make a connectivity matrix for bath spins.
- Parameters
bath (BathArray) – Array of bath spins.
r_dipole (float) – Maximum connectivity distance.
r_inner (float) – Minimum connectivity distance.
ignore (list or str, optional) – If not None, includes the names of bath spins which are ignored in the cluster generation.
max_size (int) – Maximum size of the bath before less optimal (but less memory intensive) approach is used.
- Returns
Connectivity matrix.
- Return type
crs_matrix
- connected_components(csgraph, directed=False, connection='weak', return_labels=True)
Find connected components using
scipy.sparse.csgraph
. See documentation ofscipy.sparse.csgraph.connected_components
- find_subclusters(maximum_order, graph, labels, n_components, strong=False)
Find subclusters from connectivity matrix.
- Parameters
maximum_order (int) – Maximum size of the clusters to find.
graph (csr_matrix) – Connectivity matrix.
labels (ndarray with shape (n,)) – Array of labels of the connected components.
n_components (int) – The number of connected components n.
strong (bool) – Whether to find only completely interconnected clusters (default False).
- Returns
Dictionary with keys corresponding to size of the cluster, and value corresponds to ndarray of shape (matrix, N). Here matrix is the number of clusters of given size, N is the size of the cluster. Each row contains indexes of the bath spins included in the given cluster.
- Return type
dict
- combine_clusters(cs1, cs2)
Combine two dictionaries with clusters.
- Parameters
cs1 (dict) – First cluster dictionary with keys corresponding to size of the cluster, and value corresponds to ndarray of shape (matrix, N).
cs2 (dict) – Second cluster dictionary with the same structure.
- Returns
Combined dictionary with unique clusters from both dictionaries.
- Return type
dict
- expand_clusters(sc)
Expand dict so each new cluster will include all possible additions of one more bath spin. This increases maximum size of the cluster by one.
- Parameters
sc (dict) – Initial clusters dictionary.
- Returns
Dictionary with expanded clusters.
- Return type
dict
- find_valid_subclusters(graph, maximum_order, nclusters=None, bath=None, strong=False, compute_strength=None)
Find subclusters from connectivity matrix.
- Parameters
maximum_order (int) – Maximum size of the clusters to find.
graph (csr_matrix) – Connectivity matrix.
nclusters (dict) – Dictionary which contain maximum number of clusters of the given size.
bath (BathArray) – Array of bath spins.
strong (bool) – Whether to find only completely interconnected clusters (default False).
- Returns
Dictionary with keys corresponding to size of the cluster, and value corresponds to ndarray of shape (matrix, N). Here matrix is the number of clusters of given size, N is the size of the cluster. Each row contains indexes of the bath spins included in the given cluster.
- Return type
dict
General decorators that are used to expand kernel of the RunObject
class or subclasses to the whole bath via CCE.
This module contains information about the way the cluster expansion is implemented in the package.
- cluster_expansion_decorator(_func=None, *, result_operator=<built-in function imul>, contribution_operator=<built-in function ipow>, removal_operator=<built-in function itruediv>, addition_operator=<function prod>)
Decorator for creating cluster correlation expansion of the method of
RunObject
class.- Parameters
_func (func) – Function to expand.
result_operator (func) – Operator which will combine the result of expansion (default: operator.imul).
contribution_operator (func) – Operator which will combine multiple contributions of the same cluster (default: operator.ipow) in the optimized approach.
result_operator – Operator which will combine the result of expansion (default: operator.imul).
removal_operator (func) – Operator which will remove subcluster contribution from the given cluster contribution. First argument cluster contribution, second - subcluster contribution (default: operator.itruediv).
addition_operator (func) – Group operation which will combine contributions from the different clusters into one contribution (default: np.prod).
- Returns
Expanded function.
- Return type
func
- optimized_approach(function, self, *arg, result_operator=<built-in function imul>, contribution_operator=<built-in function ipow>, **kwarg)
Optimized approach to compute cluster correlation expansion.
- Parameters
function (func) – Function to expand.
self (RunObject) – Object whose method is expanded.
*arg – list of positional arguments of the expanded function.
result_operator (func) – Operator which will combine the result of expansion (default: operator.imul).
contribution_operator (func) – Operator which will combine multiple contributions of the same cluster (default: operator.ipow).
**kwarg – Dictionary containing all keyword arguments of the expanded function.
- Returns
Expanded function.
- Return type
func
- direct_approach(function, self, *arg, result_operator=<built-in function imul>, removal_operator=<built-in function itruediv>, addition_operator=<function prod>, **kwarg)
Direct approach to compute cluster correlation expansion.
- Parameters
function (func) – Function to expand.
self (RunObject) – Object whose method is expanded.
result_operator (func) – Operator which will combine the result of expansion (default: operator.imul).
removal_operator (func) – Operator which will remove subcluster contribution from the given cluster contribution. First argument cluster contribution, second - subcluster contribution (default: operator.itruediv).
addition_operator (func) – Group operation which will combine contributions from the different clusters into one contribution (default: np.prod).
**kwarg – Dictionary containing all keyword arguments of the expanded function.
- Returns
Expanded method.
- Return type
func
- interlaced_decorator(_func=None, *, result_operator=<built-in function imul>, contribution_operator=<built-in function ipow>)
Decorator for creating interlaced cluster correlation expansion of the method of
RunObject
class.- Parameters
_func (func) – Function to expand.
result_operator (func) – Operator which will combine the result of expansion (default: operator.imul).
contribution_operator (func) – Operator which will combine multiple contributions of the same cluster (default: operator.ipow) in the optimized approach.
- Returns
Expanded method.
- Return type
func
Decorators that are used to perform bath state sampling over the kernel of RunObject
.
- generate_bath_state(bath, nbstates, seed=None, parallel=False)
Generator of the random pure \(\hat{I}_z\) bath eigenstates.
- Parameters
bath (BathArray) – Array of bath spins.
nbstates (int) – Number of random bath states to generate.
seed (int) – Optional. Seed for RNG.
parallel (bool) – True if run in parallel mode. Default False.
- Yields
List – list of the pure bath spin state vectors.
- monte_carlo_method_decorator(func)
Decorator to sample over random bath states given function.