QwaveMPS#

QwaveMPS: Package to solve waveguide QED problems using MPS

Classes#

Bins

Bin data used for analysing time-dependent quantities.

InputParams

Input / simulation parameters:

Functions#

b(→ numpy.ndarray)

Annihilation operator for observables in the truncated Fock basis.

b_dag(→ numpy.ndarray)

Creation operator for observables in the truncated Fock basis.

b_dag_l(→ numpy.ndarray)

Left creation operator for a system with two field channels in the truncated Fock basis.

b_dag_r(→ numpy.ndarray)

Right creation operator for a system with two field channels, in the truncated Fock basis.

b_l(→ numpy.ndarray)

Left annihilation operator for a system with two field channels in the truncated Fock basis.

b_pop(→ numpy.ndarray)

Single-channel photonic population operator (normalized by \(\frac{1}{\Delta t}\)).

b_pop_l(→ numpy.ndarray)

Left-channel photonic population operator (normalized by \(\frac{1}{\Delta t}\)).

b_pop_r(→ numpy.ndarray)

Right-channel photonic population operator (normalized by \(\frac{1}{\Delta t}\)).

b_r(→ numpy.ndarray)

Right annihilation operator for a system with two field channels, in the truncated Fock basis.

cite()

Print BibTeX citation for the QwaveMPS package.

correlation_2op_1t(...)

Calculates the two time correlation function \(\langle A(t_0)B(t_0+t')\rangle\) at a fixed time \(t_0\) for either single operators \(A/B\), or each operator in the lists.

correlation_2op_2t(...)

Calculates the two time correlation function \(\langle A(t)B(t+t')\rangle\) for either single operators \(A\) and \(B\), or each \(A/B\) in a_op_list/b_op_list.

correlation_4op_1t(...)

Calculates the two time correlation function \(\langle A(t_0)B(t_0+t')C(t_0+t')D(t_0)\rangle\) at a fixed time \(t_0\) for either single operators, or each operator in the lists.

correlation_4op_2t(...)

Calculates the two time correlation function \(\langle A(t)B(t+t')C(t+t')D(t)\rangle\) for either single operators \(A/B/C/D\), or each operator in the four lists.

correlation_ss_1t(→ tuple[list[numpy.ndarray], ...)

Efficient steady-state correlation calculation.

correlation_ss_2op(...)

Calculates the two time correlation function \(\langle A(t_{ss})B(t_{ss}+t')\rangle\) at a steady state value of t for either single operators, or each operator in the lists.

correlation_ss_4op(...)

Calculates the two time correlation function \(\langle A(t_{ss})B(t_{ss}+t')C(t_{ss}+t')D(t_{ss})\rangle\) at a steady state value of t for either single operators, or each operator in the lists.

correlations_1t(→ tuple[list[numpy.ndarray], ...)

General two-time correlation calculator along a single axis.

correlations_2t(→ tuple[list[numpy.ndarray], ...)

General two-time correlation calculator.

coupling(→ tuple[float, float])

Return (gamma_l, gamma_r) given a coupling specification.

delta_b(→ numpy.ndarray)

Time bin noise annihilation operator scaled by \(\sqrt{\Delta t}\) in the truncated Fock

delta_b_dag(→ numpy.ndarray)

Time bin noise creation operator scaled by \(\sqrt{\Delta t}\) in the truncated Fock

delta_b_dag_l(→ numpy.ndarray)

Left time bin noise creation operator for a system with two field channels,

delta_b_dag_r(→ numpy.ndarray)

Right time bin noise creation operator for a system with two field channels,

delta_b_l(→ numpy.ndarray)

Left time bin noise annihilation operator for a system with two field channels,

delta_b_r(→ numpy.ndarray)

Right time bin noise annihilation operator for a system with two field channels,

e(→ numpy.ndarray)

Projector onto the excited TLS state, \(|e\rangle\langle e|\).

entanglement(→ list[float])

Compute von Neumann entanglement entropy across a list of Schmidt coefficient arrays.

exp_decay_envelope(→ numpy.ndarray)

Create a exponential decay pulse envelope (unnormalized) given by the time length of the pulse

fock_pulse(→ list[numpy.ndarray])

Creates an Fock pulse input field state with a normalized pulse envelope

gaussian_envelope(→ numpy.ndarray)

Create a gaussian pulse envelope given by the time length of the pulse

hamiltonian_1tls(→ Hamiltonian)

Hamiltonian for 1 two-level system coupled to an infinite waveguide.

hamiltonian_1tls_feedback(→ Hamiltonian)

Hamiltonian for 1 two-level system in a semi-infinite waveguide with a side mirror (with feedback).

hamiltonian_2tls_mar(→ Hamiltonian)

Hamiltonian for 2 two-level systems in an infinite waveguide in the Markovian regime.

hamiltonian_2tls_nmar(→ Hamiltonian)

Hamiltonian for 2 two-level systems in an infinite waveguide in the non-Markovian regime (feedback).

input_state_generator(...)

Creates an iterator (generator) for the input field states of the waveguide.

loop_integrated_statistics(→ numpy.ndarray)

Calculates the time dependent integral of the function over all points in the feedback loop at each time point of the system evolution.

normalize_pulse_envelope(→ numpy.ndarray)

Normalizes a given pulse envelope so that the integral of the square magnitude is 1.

sigmaminus(→ numpy.ndarray)

Lowering operator for the Pauli spins, \(|g\rangle\langle e|\).

sigmaplus(→ numpy.ndarray)

Raising operator for the Pauli spins, \(|e\rangle \langle g|\).

single_time_expectation(→ numpy.ndarray)

Compute expectation values of a list of operators on a list of OC normalized bins.

spectral_intensity(→ tuple[numpy.ndarray, numpy.ndarray])

Calculate the time dependent spectral intensity from a given two time correlation function. Given a correlation function of the form \(\langle A(t)B(t+\tau)\rangle\) this computes the function

spectrum_w(→ tuple[numpy.ndarray, numpy.ndarray])

Compute the (discrete) spectrum in the long-time limit via Fourier transform

steady_state_index(→ numpy.ndarray[int])

Steady-state index helper function to find the time step

t_evol_mar(...)

Time evolution of the system without delay times (Markovian regime)

t_evol_nmar(...)

Time evolution of the system with finite delays/feedback (non-Markovian regime).

time_dependent_spectrum(→ tuple[numpy.ndarray, ...)

Calculate the time dependent spectra from a given two time correlation function. Given a correlation function of the form \(\langle A(t)B(t+\tau)\rangle\) this computes the function

tls_excited(→ numpy.ndarray)

Two level system excited state tensor.

tls_ground(→ numpy.ndarray)

Two level system ground state tensor.

tls_pop(→ numpy.ndarray)

Single TLS population operator, \(\sigma^+ \sigma^-\).

tophat_envelope(→ numpy.ndarray)

Create an unnormalized top hat pulse envelope given by the time length of the pulse.

transform_t_tau_to_t1_t2(→ numpy.ndarray)

Transforms two time correlations from a (t,tau) representation to a (t1,t2) representation.

vacuum(→ list[numpy.ndarray])

Produces an array of vacuum time bins for a given time_length.

wg_ground(→ numpy.ndarray)

Waveguide vacuum state for a single time bin.

Package Contents#

class QwaveMPS.Bins#

Bin data used for analysing time-dependent quantities.

Parameters:
  • system_states (list) – List of system bins used when calculating single time system observables.

  • output_field_states (list) – Time bins used when calculating single time field observables.

  • input_field_states (list) – List of input time bins used for calculating single time field observables incident the system resulting from the defined initial field state.

  • correlation_bins (list) – Correlation bins used when computing output field photon correlation functions at two time points.

  • schmidt (list) – Schmidt decomposition system bins usen when calculating entanglement entropy.

  • loop_field_states (list, default: None) – Tau (delay) bins used when calculating delayed field observables. This is the list of field states entering the feedback loop at each time point.

  • schmidt_tau (list, default: None) – Schmidt decomposition tau bins usen when calculating delayed entanglement entropy.

class QwaveMPS.InputParams#

Input / simulation parameters:

Parameters:
  • delta_t (float) – Time step used for time propagation.

  • tmax (float) – Maximum simulation time.

  • d_sys_total (np.ndarray) – Array describing the local physical dimensions of the system bins. Each index is associated with the size of a tensor space.

  • d_t_total (np.ndarray) – Array with dimensions for the time bins. Each index is associated with the size of a tensor space. In the case of a two directional light channel will be a list of two values.

  • bond_max (int) – Maximum MPS bond dimension (chi) to use for truncation.

  • gamma_l (float) – Coupling (decay) rates to the left and right channels respectively. (may be removed from this class in future versions)

  • gamma_r (float) – Coupling (decay) rates to the left and right channels respectively. (may be removed from this class in future versions)

  • gamma_l2 (float, default: 0) – Optional second set of coupling rates (e.g. for a second TLS). Default 0 means only one system. (may be removed from this class in future versions)

  • gamma_r2 (float, default: 0) – Optional second set of coupling rates (e.g. for a second TLS). Default 0 means only one system. (may be removed from this class in future versions)

  • tau (float, default: 0) – Delay time if modelling non-Markovian dynamics.

  • phase (float, defulat: 0) – Relative delayed phase. (may be removed from this class in future versions)

  • d_t (int) – Total size of the photonic tensor space. This is the product of the sizes of the indvidual tensor spaces.

QwaveMPS.b(params: QwaveMPS.parameters.InputParams) numpy.ndarray#

Annihilation operator for observables in the truncated Fock basis. Normalized by \(\frac{1}{\sqrt{\Delta t}}\).

Parameters:

params (InputParams) – Object that contains the simulation parameters (time step size and photonic tensor space sizes)

Returns:

oper – Annihilation operator observable.

Return type:

ndarray

QwaveMPS.b_dag(params: QwaveMPS.parameters.InputParams) numpy.ndarray#

Creation operator for observables in the truncated Fock basis. Normalized by \(\frac{1}{\sqrt{\Delta t}}\).

Parameters:

params (InputParams) – Object that contains the simulation parameters (time step size and photonic tensor space sizes)

Returns:

oper – Creation operator observable.

Return type:

ndarray

QwaveMPS.b_dag_l(params: QwaveMPS.parameters.InputParams) numpy.ndarray#

Left creation operator for a system with two field channels in the truncated Fock basis. Normalized by \(\frac{1}{\sqrt{\Delta t}}\).

Parameters:

params (InputParams) – Object that contains the simulation parameters (time step size and photonic tensor space sizes)

Returns:

oper – Creation operator observable for left channel.

Return type:

ndarray

QwaveMPS.b_dag_r(params: QwaveMPS.parameters.InputParams) numpy.ndarray#

Right creation operator for a system with two field channels, in the truncated Fock basis. Normalized by \(\frac{1}{\sqrt{\Delta t}}\).

Parameters:

params (InputParams) – Object that contains the simulation parameters (time step size and photonic tensor space sizes)

Returns:

oper – Creation operator observable for right channel.

Return type:

ndarray

QwaveMPS.b_l(params: QwaveMPS.parameters.InputParams) numpy.ndarray#

Left annihilation operator for a system with two field channels in the truncated Fock basis. Normalized by \(\frac{1}{\sqrt{\Delta t}}\).

Parameters:

params (InputParams) – Object that contains the simulation parameters (time step size and photonic tensor space sizes)

Returns:

oper – Annihilation operator observable for left channel.

Return type:

ndarray

QwaveMPS.b_pop(params: QwaveMPS.parameters.InputParams) numpy.ndarray#

Single-channel photonic population operator (normalized by \(\frac{1}{\Delta t}\)).

Parameters:

params (InputParams) – Object that contains the simulation parameters (time step size and photonic tensor space sizes)

Returns:

oper – Population operator observable.

Return type:

ndarray

QwaveMPS.b_pop_l(params: QwaveMPS.parameters.InputParams) numpy.ndarray#

Left-channel photonic population operator (normalized by \(\frac{1}{\Delta t}\)).

Parameters:

params (InputParams) – Object that contains the simulation parameters (time step size and photonic tensor space sizes)

Returns:

oper – Population operator observable for the left channel.

Return type:

ndarray

QwaveMPS.b_pop_r(params: QwaveMPS.parameters.InputParams) numpy.ndarray#

Right-channel photonic population operator (normalized by \(\frac{1}{\Delta t}\)).

Parameters:

params (InputParams) – Object that contains the simulation parameters (time step size and photonic tensor space sizes)

Returns:

oper – Population operator for the right channel

Return type:

ndarray

QwaveMPS.b_r(params: QwaveMPS.parameters.InputParams) numpy.ndarray#

Right annihilation operator for a system with two field channels, in the truncated Fock basis. Normalized by \(\frac{1}{\sqrt{\Delta t}}\).

Parameters:

params (InputParams) – Object that contains the simulation parameters (time step size and photonic tensor space sizes)

Returns:

oper – Annihlation operator observable for right channel.

Return type:

ndarray

QwaveMPS.cite()#

Print BibTeX citation for the QwaveMPS package.

QwaveMPS.correlation_2op_1t(correlation_bins: list[numpy.ndarray], a_op_list: numpy.ndarray | list[numpy.ndarray], b_op_list: numpy.ndarray | list[numpy.ndarray], t: float, params: QwaveMPS.parameters.InputParams) tuple[list[numpy.ndarray] | numpy.ndarray, numpy.ndarray]#

Calculates the two time correlation function \(\langle A(t_0)B(t_0+t')\rangle\) at a fixed time \(t_0\) for either single operators \(A/B\), or each operator in the lists. Provides list functionality as a single function call with a list of operators is much faster than individual function calls for each operator.

Parameters:
  • correlation_bins (list[ndarray]) – List of time bins with the OC in either the initial or final bin in the list.

  • a_op_list (ndarray/list[ndarray]) – Single operator, A, or a list of operators.

  • b_op_list (ndarray/list[ndarray]) – Single operator, B, or a list of operators.

  • t (float) – Fixed time point for the two time point correlation function calculation.

  • params (InputParams) – Simulation parameters

Returns:

  • correlations (list[np.ndarray]) – In the case of single operators a 1D array. In the case of a list of operators returns a list of 1D arrays, each a two time correlation function of fixed t, corresponding by index to the operators in the two operator lists. The two time correlation function is stored as f[t’], with time increments between points given by the simulation.

  • t_list (np.ndarray) – List of time points for the t’ axis.

QwaveMPS.correlation_2op_2t(correlation_bins: list[numpy.ndarray], a_op_list: numpy.ndarray | list[numpy.ndarray], b_op_list: numpy.ndarray | list[numpy.ndarray], params: QwaveMPS.parameters.InputParams, completion_print_flag: bool = True) tuple[list[numpy.ndarray] | numpy.ndarray, numpy.ndarray]#

Calculates the two time correlation function \(\langle A(t)B(t+t')\rangle\) for either single operators \(A\) and \(B\), or each \(A/B\) in a_op_list/b_op_list. Provides list functionality as a single function call with a list of operators is much faster than individual function calls for each operator.

Parameters:
  • correlation_bins (list[ndarray]) – List of time bins with the OC in either the initial or final bin in the list.

  • a_op_list (ndarray/list[ndarray]) – Single operator, A, or a list of operators.

  • b_op_list (ndarray/list[ndarray]) – Single operator, B, or a list of operators.

  • params (InputParams) – Simulation parameters

  • completion_print_flag (bool, default: True) – Prints the percent completion of the of the outer loop over t values for the calculation. Note that each loop is shorter, resulting in the percents being weighted more heavily to the start of the calculation.

Returns:

  • correlations (list[np.ndarray]) – In the case of single A and B operators a 2D array. In the case of a list of operators returns a list of 2D arrays, each a two time correlation function corresponding by index to the operators in the two operator lists. The two time correlation function is stored as f[t,t’], with non-negative t’ and time increments between points given by the simulation.

  • t_list (np.ndarray) – List of time points for the t and t’ axes.

QwaveMPS.correlation_4op_1t(correlation_bins: list[numpy.ndarray], a_op_list: numpy.ndarray | list[numpy.ndarray], b_op_list: numpy.ndarray | list[numpy.ndarray], c_op_list: numpy.ndarray | list[numpy.ndarray], d_op_list: numpy.ndarray | list[numpy.ndarray], t: float, params: QwaveMPS.parameters.InputParams) tuple[list[numpy.ndarray] | numpy.ndarray, numpy.ndarray]#

Calculates the two time correlation function \(\langle A(t_0)B(t_0+t')C(t_0+t')D(t_0)\rangle\) at a fixed time \(t_0\) for either single operators, or each operator in the lists. Provides list functionality as a single function call with a list of operators is much faster than individual function calls for each operator.

Parameters:
  • correlation_bins (list[ndarray]) – List of time bins with the OC in either the initial or final bin in the list.

  • a_op_list (ndarray/list[ndarray]) – Single operator, A, or a list of operators.

  • b_op_list (ndarray/list[ndarray]) – Single operator, B, or a list of operators.

  • c_op_list (ndarray/list[ndarray]) – Single operator, C, or a list of operators.

  • d_op_list (ndarray/list[ndarray]) – Single operator, D, or a list of operators.

  • t (float) – Fixed time point for the two time point correlation function calculation.

  • params (InputParams) – Simulation parameters

Returns:

  • correlations (list[np.ndarray]) – In the case of single operators a 1D array. In the case of a list of operators returns a list of 1D arrays, each a two time correlation function of fixed t, corresponding by index to the operators in the two operator lists. The two time correlation function is stored as f[t’], with time increments between points given by the simulation.

  • t_list (np.ndarray) – List of time points for the t’ axis.

QwaveMPS.correlation_4op_2t(correlation_bins: list[numpy.ndarray], a_op_list: numpy.ndarray | list[numpy.ndarray], b_op_list: numpy.ndarray | list[numpy.ndarray], c_op_list: numpy.ndarray | list[numpy.ndarray], d_op_list: numpy.ndarray | list[numpy.ndarray], params: QwaveMPS.parameters.InputParams, completion_print_flag: bool = True) tuple[list[numpy.ndarray] | numpy.ndarray, numpy.ndarray]#

Calculates the two time correlation function \(\langle A(t)B(t+t')C(t+t')D(t)\rangle\) for either single operators \(A/B/C/D\), or each operator in the four lists. Provides list functionality as a single function call with a list of operators is much faster than individual function calls for each operator.

Parameters:
  • correlation_bins (list[ndarray]) – List of time bins with the OC in either the initial or final bin in the list.

  • a_op_list (ndarray/list[ndarray]) – Single operator, A, or a list of operators.

  • b_op_list (ndarray/list[ndarray]) – Single operator, B, or a list of operators.

  • c_op_list (ndarray/list[ndarray]) – Single operator, C, or a list of operators.

  • d_op_list (ndarray/list[ndarray]) – Single operator, D, or a list of operators.

  • params (InputParams) – Simulation parameters

  • completion_print_flag (bool, default: True) – Prints the percent completion of the of the outer loop over t values for the calculation. Note that each loop is shorter, resulting in the percents being weighted more heavily to the start of the calculation.

Returns:

  • correlations (list[np.ndarray]) – In the case of single operators a 2D array. In the case of a list of operators returns a list of 2D arrays, each a two time correlation function corresponding by index to the operators in the two operator lists. The two time correlation function is stored as f[t,t’], with non-negative t’ and time increments between points given by the simulation.

  • t_list (np.ndarray) – List of time points for the t and t’ axes.

QwaveMPS.correlation_ss_1t(correlation_bins: list[numpy.ndarray], output_field_states: list[numpy.ndarray], ops_same_time: list[numpy.ndarray], ops_two_time: list[numpy.ndarray], params: QwaveMPS.parameters.InputParams, tol: float = 1e-05, window: int = 10, t_steady: float = None) tuple[list[numpy.ndarray], numpy.ndarray, float]#

Efficient steady-state correlation calculation. This computes time differences starting from a convergence index (steady-state index). It returns a list of the 1D correlation arrays corresponding to the operator list, a list of tau points, and the initial t point at which steady state is considered.

Parameters:
  • correlation_bins (list[np.ndarray]) – Correlation bins of the outputfield states used to determine multi-timepoint correlation functions of the output field.

  • output_field_states (list[np.ndarray]) – OC normalized output field states.

  • ops_same_time ([ndarray]) – List of operators of which correlation functions should be calculated in the case that tau=0 (same time). These should exist in a single time-bin tensor space.

  • ops_two_time ([ndarray]) – List of operators of which correlation functions should be calculated in the case that tau > 0. These should be ordered in a corresponding order to ops_same_time and should exist in a tensor space that is the outer product of two time bin tensor spaces, with the right space corresponding to the greater time.

  • params (InputParams) – Simulation parameters

  • window (int, default: 10) – Number of recent points to analyze when determining the steady state time.

  • tol (float, default: 1e-5) – Maximum deviation allowed in the final window for the steady state time.

  • t_steady (float, default: None) – User defined steady state time. If not provided, steady state is determined by convergence of the same time expectation values of the observables.

Returns:

  • correlations (list[ndarray]) – A list of 1D arrays, each a two time correlation function of fixed t at steady state, corresponding by index to the operators in the two operator lists. The two time correlation function is stored as f[t’], with time increments between points given by the simulation.

  • t_list (ndarray) – List of time points for the t’ axis.

  • t_ss (float) – Time that steady state is reached.

QwaveMPS.correlation_ss_2op(correlation_bins: list[numpy.ndarray], output_field_states: list[numpy.ndarray], a_op_list: numpy.ndarray | list[numpy.ndarray], b_op_list: numpy.ndarray | list[numpy.ndarray], params: QwaveMPS.parameters.InputParams, tol: float = 1e-05, window: int = 20, t_steady: float = None) tuple[list[numpy.ndarray] | numpy.ndarray, numpy.ndarray, float]#

Calculates the two time correlation function \(\langle A(t_{ss})B(t_{ss}+t')\rangle\) at a steady state value of t for either single operators, or each operator in the lists. Provides list functionality as a single function call with a list of operators is much faster than individual function calls for each operator. In that case calculates the steady states correlation from the greatest steady state time of the operators.

Parameters:
  • correlation_bins (list[ndarray]) – List of time bins with the OC in either the initial or final bin in the list.

  • a_op_list (ndarray/list[ndarray]) – Single operator, A, or a list of operators.

  • b_op_list (ndarray/list[ndarray]) – Single operator, B, or a list of operators.

  • params (InputParams) – Simulation parameters

  • tol (float, default: 1e-5) – The tolerance for which convergence of the operators is determined. Used to find the steady state time.

  • window (int, default: 20) – Number of recent points to analyze when determining the steady state time.

  • t_steady (float, default: None) – User defined steady state time. If not provided, steady state is determined by convergence of the same time expectation values of the observables.

Returns:

  • correlations (list[ndarray]) – In the case of single operators a 1D array. In the case of a list of operators returns a list of 1D arrays, each a two time correlation function of fixed t at steady state, corresponding by index to the operators in the two operator lists. The two time correlation function is stored as f[t’], with time increments between points given by the simulation.

  • t_list (ndarray) – List of time points for the t’ axis.

  • t_ss (float) – Time that steady state is reached.

QwaveMPS.correlation_ss_4op(correlation_bins: list[numpy.ndarray], output_field_states: list[numpy.ndarray], a_op_list: numpy.ndarray | list[numpy.ndarray], b_op_list: numpy.ndarray | list[numpy.ndarray], c_op_list: numpy.ndarray | list[numpy.ndarray], d_op_list: numpy.ndarray | list[numpy.ndarray], params: QwaveMPS.parameters.InputParams, tol: float = 1e-05, window: int = 20, t_steady: float = None) tuple[list[numpy.ndarray] | numpy.ndarray, numpy.ndarray, float]#

Calculates the two time correlation function \(\langle A(t_{ss})B(t_{ss}+t')C(t_{ss}+t')D(t_{ss})\rangle\) at a steady state value of t for either single operators, or each operator in the lists. Provides list functionality as a single function call with a list of operators is much faster than individual function calls for each operator. In that case calculates the steady states correlation from the greatest steady state time of the operators.

Parameters:
  • correlation_bins (list[ndarray]) – List of time bins with the OC in either the initial or final bin in the list.

  • a_op_list (ndarray/list[ndarray]) – Single operator, A, or a list of operators.

  • b_op_list (ndarray/list[ndarray]) – Single operator, B, or a list of operators.

  • c_op_list (ndarray/list[ndarray]) – Single operator, C, or a list of operators.

  • d_op_list (ndarray/list[ndarray]) – Single operator, D, or a list of operators.

  • params (InputParams) – Simulation parameters

  • tol (float, default: 1e-5) – The tolerance for which convergence of the operators is determined. Used to find the steady state time.

  • window (int, default: 20) – Number of recent points to analyze when determining the steady state time.

  • t_steady (float, default: None) – User defined steady state time. If not provided, steady state is determined by convergence of the same time expectation values of the observables.

Returns:

  • correlations (list[ndarray]) – In the case of single operators a 1D array. In the case of a list of operators returns a list of 1D arrays, each a two time correlation function of fixed t at steady state, corresponding by index to the operators in the two operator lists. The two time correlation function is stored as f[t’], with time increments between points given by the simulation.

  • t_list (ndarray) – List of time points for the t’ axis.

  • t_ss (float) – Time that steady state is reached.

QwaveMPS.correlations_1t(correlation_bins: list[numpy.ndarray], ops_same_time: list[numpy.ndarray], ops_two_time: list[numpy.ndarray], t: float, params: QwaveMPS.parameters.InputParams) tuple[list[numpy.ndarray], numpy.ndarray]#

General two-time correlation calculator along a single axis. Take in list of time ordered normalized (with OC) time bins at position of relevance. Calculate a list of arbitrary two time point correlation functions at t and t+tau for nonnegative t’.

Parameters:
  • time_bin_list (list[ndarray]) – List of time bins with the OC in either the initial or final bin in the list.

  • ops_same_time (list[ndarray]) – List of operators of which correlation functions should be calculated in the case that t’=0 (same time). These should exist in a single time-bin tensor space.

  • ops_two_time (list[ndarray]) – List of operators of which correlation functions should be calculated in the case that t’ > 0. These should be ordered in a corresponding order to ops_same_time and should exist in a tensor space that is the outer product of two time bin tensor spaces, with the right space corresponding to the greater time.

  • t (float) – Time point for fixed t at which to take the two time point correlation.

  • params (InputParams) – Simulation parameters

Returns:

  • correlations (ndarray[ndarray[complex]]) – List of 1D arrays, each a two time correlation function corresponding by index to the operators in ops_same_time and ops_two_time. The two time correlation function is stored as f[t,t’], with non-negative t’ and time increments between points given by the simulation.

  • ts_correlation (ndarray[float]) – List of time points for the t’ axis at which the two time point correlation functions are taken.

QwaveMPS.correlations_2t(correlation_bins: list[numpy.ndarray], ops_same_time: list[numpy.ndarray], ops_two_time: list[numpy.ndarray], params: QwaveMPS.parameters.InputParams, completion_print_flag: bool = False) tuple[list[numpy.ndarray], numpy.ndarray]#

General two-time correlation calculator. Take in list of time ordered normalized (with OC) time bins at position of relevance. Calculate a list of arbitrary two time point correlation functions at t and t+t’ for nonnegative t’.

Parameters:
  • time_bin_list ([ndarray]) – List of time bins with the OC in either the initial or final bin in the list.

  • ops_same_time ([ndarray]) – List of operators of which correlation functions should be calculated in the case that t’=0 (same time). These should exist in a single time-bin tensor space.

  • ops_two_time ([ndarray]) – List of operators of which correlation functions should be calculated in the case that t’ > 0. These should be ordered in a corresponding order to ops_same_time and should exist in a tensor space that is the outer product of two time bin tensor spaces, with the right space corresponding to the greater time.

  • params (InputParams) – Simulation parameters

  • completion_print_flag (bool, default=True) – Flag to print completion loop number percent of the calculation (note this is not the percent completion, and later loops complete faster than earlier ones).

Returns:

  • result (list[np.ndarray]) – List of 2D arrays, each a two time correlation function corresponding by index to the operators in ops_same_time and ops_two_time. The two time correlation function is stored as f[t,t’], with non-negative t’ and time increments between points given by the simulation.

  • correaltion_times (np.ndarray[float]) – List of time points for the t and t’ axes for the calculated correlation functions.

QwaveMPS.coupling(coupl: str = 'symmetrical', gamma: float = 1, gamma_r=None, gamma_l=None) tuple[float, float]#

Return (gamma_l, gamma_r) given a coupling specification.

It can be ‘symmetrical’, ‘chiral_r’, ‘chiral_l’, ‘other’ For ‘other’, provide gamma_l and gamma_r explicitly.

Parameters:
  • coupl ({'symmetrical', 'chiral_r', 'chiral_l', 'other'}, default: 'symmetrical') – Coupling option.

  • gamma (float, default:1) – Total coupling. Code in units of coupling, hence, the default is 1.

  • gamma_r (None/float, default: None) – Left coupling. If coupl = ‘other’ define explicitly.

  • gamma_l (None/float, default: None) – Right coupling. If coupl = ‘other’ define explicitly.

Returns:

gamma_l,gamma_r – Values of the left and right coupling

Return type:

tuple[float,float]

QwaveMPS.delta_b(delta_t: float, d_t: int = 2) numpy.ndarray#

Time bin noise annihilation operator scaled by \(\sqrt{\Delta t}\) in the truncated Fock basis.

Parameters:
  • delta_t (float) – Time step for system evolution.

  • d_t (int, default: 2) – Size of the truncated field Hilbert space

Returns:

oper – Time bin noise creation operator.

Return type:

ndarray

QwaveMPS.delta_b_dag(delta_t: float, d_t: int = 2) numpy.ndarray#

Time bin noise creation operator scaled by \(\sqrt{\Delta t}\) in the truncated Fock basis.

Parameters:
  • delta_t (float) – Time step for system evolution.

  • d_t (int, default: 2) – Size of the truncated field Hilbert space

Returns:

oper – Time bin noise creation operator.

Return type:

ndarray

QwaveMPS.delta_b_dag_l(delta_t: float, d_t_total: numpy.ndarray) numpy.ndarray#

Left time bin noise creation operator for a system with two field channels, scaled by \(\sqrt{\Delta t}\) in the truncated Fock basis.

Parameters:
  • delta_t (float) – Time step for system evolution.

  • d_t_total (ndarray) – List of sizes of the photonic Hilbert spaces (left and right channels).

Returns:

oper – Left time bin noise creation operator.

Return type:

ndarray

QwaveMPS.delta_b_dag_r(delta_t: float, d_t_total: numpy.ndarray) numpy.ndarray#

Right time bin noise creation operator for a system with two field channels, scaled by \(\sqrt{\Delta t}\) in the truncated Fock basis.

Parameters:
  • delta_t (float) – Time step for system evolution.

  • d_t_total (ndarray) – List of sizes of the photonic Hilbert spaces (left and right channels).

Returns:

oper – Right time bin noise creation operator.

Return type:

ndarray

QwaveMPS.delta_b_l(delta_t: float, d_t_total: numpy.ndarray) numpy.ndarray#

Left time bin noise annihilation operator for a system with two field channels, scaled by \(\sqrt{\Delta t}\) in the truncated Fock basis.

Parameters:
  • delta_t (float) – Time step for system evolution.

  • d_t_total (ndarray) – List of sizes of the photonic Hilbert spaces (left and right channels).

Returns:

oper – Left time bin noise annihilation operator.

Return type:

ndarray

QwaveMPS.delta_b_r(delta_t: float, d_t_total: numpy.ndarray) numpy.ndarray#

Right time bin noise annihilation operator for a system with two field channels, scaled by \(\sqrt{\Delta t}\) in the truncated Fock basis.

Parameters:
  • delta_t (float) – Time step for system evolution.

  • d_t_total (ndarray) – List of sizes of the photonic Hilbert spaces.

Returns:

oper – Right time bin noise annihilation operator (left and right channels).

Return type:

ndarray

QwaveMPS.e(d_sys: int = 2) numpy.ndarray#

Projector onto the excited TLS state, \(|e\rangle\langle e|\).

Parameters:

d_sys (int, default: 2) – Size of the Hilbert space of the matter system. For a two level system have d_sys=2.

Returns:

oper – The excited state projector

Return type:

ndarray

QwaveMPS.entanglement(sch: list[numpy.ndarray]) list[float]#

Compute von Neumann entanglement entropy across a list of Schmidt coefficient arrays.

Parameters:

sch (list[np.ndarray]) – List of Schmidt coefficient arrays (s) for each bipartition.

Returns:

time_dependent_entanglement – Entanglement entropies computed as \(-\sum(p\log_2 p)\) where \(p = s^2\).

Return type:

list[float]

QwaveMPS.exp_decay_envelope(pulse_time: float, params: QwaveMPS.parameters.InputParams, decay_rate: float, decay_center: float = 0) numpy.ndarray#

Create a exponential decay pulse envelope (unnormalized) given by the time length of the pulse and the decay rate and decay center parameters.

Parameters:
  • pulse_time (float) – Duration time of the pulse (units of inverse coupling).

  • params (InputParams) – Class containing the input parameters

  • decay_rate (float) – Decay rate of the exponential.

  • decay_center (float) – Time center/offset of the exponential decay function.

Returns:

pulse_envelope – List of amplitude values of the pulse envelope.

Return type:

list[float]

QwaveMPS.fock_pulse(pulse_env: list[float], pulse_time: float, photon_num: int, params: QwaveMPS.parameters.InputParams, direction: str = 'R', bond0: int = 1) list[numpy.ndarray]#

Creates an Fock pulse input field state with a normalized pulse envelope

Parameters:
  • pulse_env (list[float]) – Time dependent pulse envelope for the incident pulse (can be unnormalized). If None, uses a tophat pulse for the duration of the pulse_time.

  • pulse_time (float) – Time length of the pulse (units of inverse coupling). If the pulse envelope is of greater length it will be truncated from the tail.

  • photon_num (int) – Incident photon number.

  • params (InputParams) – Class containing the input parameters

  • direction ({'L','R'}, default: 'R') – Flag to dictate the direction of the propagating pulse. Ignored if only a single photonic channel (chiral) is present.

  • bond0 (int, default: 1) – Default bond dimension of bins.

Returns:

fock_pulse – A list of the incident time bins of the Fock pulse, with the first bin in index 0. Further uncorrelated fields can be appended to the end of the list.

Return type:

list[ndarray]

Examples

QwaveMPS.gaussian_envelope(pulse_time: float, params: QwaveMPS.parameters.InputParams, gaussian_width: float, gaussian_center: float) numpy.ndarray#

Create a gaussian pulse envelope given by the time length of the pulse and the mean and standard deviation parameters.

Parameters:
  • pulse_time (float) – Duration time of the pulse (units of inverse coupling).

  • params (InputParams) – Class containing the input parameters

  • gaussian_width (float) – Variance of the gaussian (units of inverse coupling).

  • gaussian_center (float) – Mean of the gaussian (units of inverse coupling).

Returns:

pulse_envelope – List of amplitude values of the pulse envelope.

Return type:

np.ndarray[float]

QwaveMPS.hamiltonian_1tls(params: QwaveMPS.parameters.InputParams, omega: float | QwaveMPS.operators.np.ndarray = 0, delta: float = 0) Hamiltonian#

Hamiltonian for 1 two-level system coupled to an infinite waveguide.

The returned Hamiltonian includes:
  • A classical pump term (omega) acting on the TLS, \(\Omega(\sigma^+ + \sigma^-)\)

  • A detuning term for the TLS, \(\delta |e\rangle\langle e|\).

  • Interaction terms between the TLS and left/right photonic modes.

Parameters:
  • params (InputParams) – Class containing the input parameters.

  • omega (float/np.ndarray, default: 0) – Classical pump amplitude. If a float is provided (CW pump) a single Hamiltonian ndarray is returned. If a 1D np.ndarray is given (pulsed light), the function returns a callable hm_total(t_k) that yields the Hamiltonian at discrete time index t_k using omega[t_k].

  • delta (float, optional) – Detuning between the pump and two-level system transition frequency.

Returns:

Hamiltonian – Hamiltonian as a numpy.ndarray (time-independent drive) or a callable that accepts a time index and returns the Hamiltonian (time-dependent drive).

Return type:

np.ndarray | Callable[[int], np.ndarray]

QwaveMPS.hamiltonian_1tls_feedback(params: QwaveMPS.parameters.InputParams, omega: float | QwaveMPS.operators.np.ndarray = 0, delta: float = 0) Hamiltonian#

Hamiltonian for 1 two-level system in a semi-infinite waveguide with a side mirror (with feedback).

The returned Hamiltonian includes:
  • A classical pump term (omega) acting on the TLS, \(\Omega(\sigma^+ + \sigma^-)\)

  • A detuning term for the TLS, \(\delta |e\rangle\langle e|\).

  • Interaction terms between the TLS and a single photonic mode (on the present and feedback bins).

Parameters:
  • params (InputParams) – Class containing the input parameters. (It must include the phase)

  • omega (float or np.ndarray, default: 0) – Classical pump amplitude. If a float is provided (CW pump) a single Hamiltonian ndarray is returned. If a 1D np.ndarray is given (pulsed light), the function returns a callable hm_total(t_k) that yields the Hamiltonian at discrete time index t_k using omega[t_k].

  • delta (float, default: 0) – Detuning between the pump and TLS transition frequency.

Returns:

Hamiltonian – Hamiltonian as a numpy.ndarray (time-independent drive) or a callable that accepts a time index and returns the Hamiltonian (time-dependent drive).

Return type:

np.ndarray | Callable[[int], np.ndarray]

QwaveMPS.hamiltonian_2tls_mar(params: QwaveMPS.parameters.InputParams, omega1: float | QwaveMPS.operators.np.ndarray = 0, delta1: float = 0, omega2: float | QwaveMPS.operators.np.ndarray = 0, delta2: float = 0) Hamiltonian#

Hamiltonian for 2 two-level systems in an infinite waveguide in the Markovian regime.

The returned Hamiltonian includes:
  • Classical pump terms (omega1/omega2) acting on TLS1/TLS2, \(\Omega_i(\sigma_i^+ + \sigma_i^-)\)

  • A detuning term delta1/delta2 for TLS1/TLS2, \(\delta_i |e\rangle_i\langle e|_i\).

  • Interaction terms between the TLSs and left/right photonic modes.

Parameters:
  • params (InputParams) – Class containing the input parameters.

  • omega1 (float/np.ndarray, default: 0) – Drive for two-level system 1 (can be a float for CW pumps or a time-dependent array for pulsed light).

  • delta1 (float, default: 0) – Detuning for two-level system 1.

  • omega2 (float/np.ndarray, default: 0) – Drive for two-level system 2 (can be a float for CW pumps or a time-dependent array for pulsed light).

  • delta2 (float, default: 0) – Detuning for two-level system 1.

Returns:

Hamiltonian – Hamiltonian as a numpy.ndarray (time-independent drive) or a callable that accepts a time index and returns the Hamiltonian (time-dependent drive).

Return type:

np.ndarray | Callable[[int], np.ndarray]

QwaveMPS.hamiltonian_2tls_nmar(params: QwaveMPS.parameters.InputParams, omega1: float | QwaveMPS.operators.np.ndarray = 0, delta1: float = 0, omega2: float | QwaveMPS.operators.np.ndarray = 0, delta2: float = 0) Hamiltonian#

Hamiltonian for 2 two-level systems in an infinite waveguide in the non-Markovian regime (feedback).

The returned Hamiltonian includes:
  • Classical pump terms (omega1/omega2) acting on TLS1/TLS2, \(\Omega_i(\sigma_i^+ + \sigma_i^-)\)

  • A detuning term delta1/delta2 for TLS1/TLS2, \(\delta_i |e\rangle_i\langle e|_i\).

  • Interaction terms between the two-level systems and left/right photonic modes (on the present and feedback bins).

Parameters:
  • params (InputParams) – Class containing the input parameters.

  • omega1 (float/np.ndarray, default: 0) – Drive for two-level system 1 (can be a float for CW pumps or a time-dependent array for pulsed light).

  • delta1 (float, default: 0) – Detuning for two-level system 1.

  • omega2 (float/np.ndarray, default: 0) – Drive for two-level system 2 (can be a float for CW pumps or a time-dependent array for pulsed light).

  • delta2 (float, default: 0) – Detuning for two-level system 1.

Returns:

Hamiltonian – Hamiltonian as a numpy.ndarray (time-independent drive) or a callable that accepts a time index and returns the Hamiltonian (time-dependent drive).

Return type:

np.ndarray | Callable[[int], np.ndarray]

QwaveMPS.input_state_generator(d_t_total: list[int], input_bins: list[numpy.ndarray] = None, bond0: int = 1, default_state=None) collections.abc.Iterator[numpy.ndarray]#

Creates an iterator (generator) for the input field states of the waveguide.

Parameters:
  • d_t_total (list[int]) – List of sizes of the photonic Hilbert spaces.

  • input_bins (list[np.ndarray], default: None) – List of time bins describing the input field state.

  • bond0 (int, default: 1) – Size of the initial bond dimension.

  • default_state (ndarray, default: None) – Default time bin state yielded as an input state after all input_bins are exhusted. If None then vacuum states are yielded.

Returns:

gen – A generator for the input field time bins.

Return type:

Generator

QwaveMPS.loop_integrated_statistics(time_dependent_func: numpy.ndarray[complex], params: QwaveMPS.parameters.InputParams) numpy.ndarray#

Calculates the time dependent integral of the function over all points in the feedback loop at each time point of the system evolution. This is a moving windowed integral over the function, with the window size of length tau (the feedback time) For t<tau assumes time points yet to be reached in the loop will initially be 0.

Parameters:
  • time_dependent_func (np.ndarray[complex]) – List of values for the time dependent function to be integrated over the loop at each time point.

  • params (InputParams) – Simulation parameters

Returns:

observable_integrated_in_loop – List of values for the integration of time_dependent_func over the feedback loop at each time point. This is a moving integral over the function, for a window of length tau.

Return type:

np.ndarray

QwaveMPS.normalize_pulse_envelope(delta_t: float, pulse_env: numpy.ndarray) numpy.ndarray#

Normalizes a given pulse envelope so that the integral of the square magnitude is 1.

Parameters:
  • delta_t (float) – Time step size for the simulation.

  • pulse_env (np.ndarray[float]) – Time dependent pulse envelope that is being normalized.

Returns:

pulse_env – The normalized time dependent pulse envelope.

Return type:

np.ndarray[float]

QwaveMPS.sigmaminus() numpy.ndarray#

Lowering operator for the Pauli spins, \(|g\rangle\langle e|\).

Returns:

sigma_minus – The Pauli spin lowering operator.

Return type:

ndarray

QwaveMPS.sigmaplus() numpy.ndarray#

Raising operator for the Pauli spins, \(|e\rangle \langle g|\).

Returns:

sigma_plus – The Pauli spin raising operator.

Return type:

ndarray

QwaveMPS.single_time_expectation(normalized_bins: list[numpy.ndarray], ops_list: numpy.ndarray | list[numpy.ndarray]) numpy.ndarray#

Compute expectation values of a list of operators on a list of OC normalized bins.

Parameters:
  • normalized_bins (list[ndarray]) – List of OC normalized bins in order of time to have localized expectation values taken.

  • ops_list (ndarray/list[ndarray]) – List of operators or a single operator to take expectation values. Each operator must be compatible with the bin physical space.

Returns:

expectation_values – In the case of a single operator returns a list of expectation values for each time point. For a list of operators returns a 2D array shaped (len(ops_list), len(normalized_bins)) with expectation values for each operator at each time.

Return type:

np.ndarray[complex] | np.ndarray[np.ndarray[complex]]

QwaveMPS.spectral_intensity(correlation_matrix: numpy.ndarray, input_params: QwaveMPS.parameters.InputParams, padding: int = 0, hanning_filter: bool = False, taper_length: int = 16) tuple[numpy.ndarray, numpy.ndarray]#

Calculate the time dependent spectral intensity from a given two time correlation function. Given a correlation function of the form \(\langle A(t)B(t+\tau)\rangle\) this computes the function

\[I(\omega, t) = \int_0^\infty d\tau \langle A(t)B(t+\tau) \rangle e^{i\Delta_{\omega p}\tau}\]
Parameters:
  • correlation_matrix (np.ndarray) – Computed two time correlation matrix used for the calculation of the spectral intensity.

  • input_params (InputParams) – Input parameters of the simulation

  • padding (int, default=0) – Number of 0’s added to the Fourier transform as padding for smoother results.

  • hanning_filter (bool, default=False) – Determines whether or not a Hanning filter is used to smooth the decay at the end of the function for a smoother result.

  • taper_length (int, default=16) – Determines the number of time points from the end of the data on which the Hanning filter is applied. Only relevant if hanning_filter is True.

Returns:

  • spectral_intensity (np.ndarray) – The computed time dependent spectral intensity of the given correlation function.

  • w_list (np.ndarray) – List of frequencies associated with the calculated spectral intensity.

QwaveMPS.spectrum_w(delta_t: float, g1_list: numpy.ndarray) tuple[numpy.ndarray, numpy.ndarray]#

Compute the (discrete) spectrum in the long-time limit via Fourier transform of the two-time first-order correlation (steady-state solution).

Parameters:
  • delta_t (float) – Time step used in the simulation; used to set frequency sampling.

  • g1_list (np.ndarray) – Steady-state first order correlation.

Returns:

  • s_w (np.ndarray) – Spectrum in the long-time limit (steady state solution)

  • wlist (np.ndarray) – Corresponding frequency list.

QwaveMPS.steady_state_index(output_field_states: list[numpy.ndarray], tol: float = 1e-05, window: int = 10) numpy.ndarray[int]#

Steady-state index helper function to find the time step when the steady state is reached in the single time dynamics of each operator.

Parameters:
  • output_field_states (list[np.ndarray]) – List of OC normalized output_field_states/bins.

  • tol (float, default: 1e-5) – Maximum deviation allowed in the final window

  • window (int, default: 10) – Number of recent points to analyze

Returns:

steady_state_index – The index of the start of the steady window for the output field.

Return type:

int

QwaveMPS.t_evol_mar(ham: QwaveMPS.hamiltonians.Hamiltonian, i_s0: QwaveMPS.operators.np.ndarray, i_n0: QwaveMPS.operators.np.ndarray, params: QwaveMPS.operators.InputParams) tuple[list[QwaveMPS.operators.np.ndarray], list[QwaveMPS.operators.np.ndarray]]#

Time evolution of the system without delay times (Markovian regime)

Parameters:
  • ham (ndarray or callable) – Either a fixed evolution operator/tensor or a callable returning the evolution operator for time-step k: ham(k).

  • i_s0 (ndarray) – Initial system bin (tensor).

  • i_n0 (ndarray) – Initial field bin. Seed for the input time-bin generator.

  • params (InputParams) – Class containing the input parameters (contains delta_t, tmax, bond, d_t_total, d_sys_total).

Returns:

results

containing:
  • sys_b: list of system bins

  • time_b: list of time bins

  • cor_b: list of tensors used for correlations

  • schmidt: list of Schmidt coefficient arrays (for entanglement calculation)

Return type:

Bins (from parameters.py)

QwaveMPS.t_evol_nmar(ham: QwaveMPS.hamiltonians.Hamiltonian, i_s0: QwaveMPS.operators.np.ndarray, i_n0: QwaveMPS.operators.np.ndarray, params: QwaveMPS.operators.InputParams) tuple[list[QwaveMPS.operators.np.ndarray], list[QwaveMPS.operators.np.ndarray], list[QwaveMPS.operators.np.ndarray]]#

Time evolution of the system with finite delays/feedback (non-Markovian regime). Requires tau to be at least delta_t.

Parameters:

ham (ndarray/callable) –

Either a fixed evolution operator/tensor or a callable returning the

evolution operator for time-step k: ham(k).

i_s0ndarray

Initial system bin (tensor).

i_n0: ndarray

Initial field bin. Seed for the input time-bin generator.

paramsInputParams

Class containing the input parameters (contains delta_t, tmax, bond, d_t_total, d_sys_total, tau.).

Returns:

Bins

containing:
  • sys_b: list of system bins

  • time_b: list of time bins

  • tau_b: list of feedback bins

  • cor_b: list of tensors used for correlations

  • schmidt, schmidt_tau: lists of Schmidt coefficient arrays

Return type:

Dataclass (from parameters.py)

QwaveMPS.time_dependent_spectrum(correlation_matrix: numpy.ndarray, input_params: QwaveMPS.parameters.InputParams, w_list: numpy.ndarray = None, padding: int = 0) tuple[numpy.ndarray, numpy.ndarray]#

Calculate the time dependent spectra from a given two time correlation function. Given a correlation function of the form \(\langle A(t)B(t+\tau)\rangle\) this computes the function

\[S(\omega, t) = \int_0^t dt' \int_0^{t-t'} d\tau \langle A(t)B(t+\tau) \rangle e^{i\Delta_{\omega p}\tau}\]
Parameters:
  • correlation_matrix (np.ndarray) – Computed two time correlation matrix used for the calculation of the spectral intensity.

  • input_params (InputParams) – Input parameters of the simulation

  • w_list (np.ndarray, default=None) – Frequency points at which the time dependent spectrum is calculated. If None, generates the frequency list using np.fft.fftfreq() on the length of the correlation_matrix.

  • padding (int, default: 0) – Padding added to the frequency domain.

Returns:

  • spectrum (np.ndarray) – The computed time dependent spectrum of the given correlation function.

  • w_list (np.ndarray) – List of frequencies associated with the calculated time dependent spectrum.

QwaveMPS.tls_excited(bond0: int = 1) numpy.ndarray#

Two level system excited state tensor.

Parameters:

bond0 (int, default: 1) – Initial size of the bond dimension.

Returns:

state – Excited state of the two level system.

Return type:

ndarray

QwaveMPS.tls_ground(bond0: int = 1) numpy.ndarray#

Two level system ground state tensor.

Parameters:

bond0 (int, default: 1) – Initial size of the bond dimension.

Returns:

state – Ground state of the two level system.

Return type:

ndarray

QwaveMPS.tls_pop(d_sys: int = 2) numpy.ndarray#

Single TLS population operator, \(\sigma^+ \sigma^-\).

Parameters:

d_sys (int, default: 2) – Size of the Hilbert space of the matter system. For a two level system have d_sys=2.

Returns:

pop_operator – Population operator for a TLS

Return type:

ndarray

QwaveMPS.tophat_envelope(pulse_time: float, params: QwaveMPS.parameters.InputParams) numpy.ndarray#

Create an unnormalized top hat pulse envelope given by the time length of the pulse.

Parameters:
  • pulse_time (float) – Duration time of the pulse (units of inverse coupling).

  • params (InputParams) – Class containing the input parameters.

Returns:

pulse_envelope – List of amplitude values of the pulse envelope.

Return type:

list[float]

QwaveMPS.transform_t_tau_to_t1_t2(positive_tau_results: numpy.ndarray, negative_tau_results: numpy.ndarray = None) numpy.ndarray#

Transforms two time correlations from a (t,tau) representation to a (t1,t2) representation. Takes the computed correlation function with operators ordered for the positive and negative tau (t>t+tau) cases. Note that this truncates the result to maintain the same overall shape by having t2 have the same domain as t1 (truncates cases where t+tau>t1_max). If only given one matrix assumes symmetry over the tau axis.

Parameters:
  • positive_tau_results (np.ndarray) – Computed two time correlation function in the case of operators ordered for positive tau data.

  • negative_tau_results (np.ndarray, default: None) – Computed two time correlation function in the case of operators ordered for negative tau data. If None, uses the positive_tau_results, treating the observable as symmetric over the tau axis.

Returns:

transformed_t1_t2_data – Truncated data with (t1,t2) axes.

Return type:

np.ndarray

QwaveMPS.vacuum(time_length: float, params: QwaveMPS.parameters.InputParams) list[numpy.ndarray]#

Produces an array of vacuum time bins for a given time_length.

Parameters:
  • time_length (float) – Length of the vacuum pulse (units of inverse coupling).

  • params (InputParams) – Class containing the input parameters.

Returns:

state – List of vacuum states for time_length.

Return type:

list[np.ndarray]

QwaveMPS.wg_ground(d_t: int, bond0: int = 1) numpy.ndarray#

Waveguide vacuum state for a single time bin.

Parameters:
  • d_t (int) – Size of the truncated Hilbert space of the light field.

  • bond0 (int, default: 1) – Initial size of the bond dimension.

Returns:

state – Waveguide vacuum state.

Return type:

ndarray