Environments

Bosonic Environments

class BosonicEnvironment(T: float = None, tag: Any = None)[source]

The bosonic environment of an open quantum system. It is characterized by its spectral density and temperature or, equivalently, its power spectrum or its two-time auto-correlation function.

Use one of the classmethods from_spectral_density, from_power_spectrum or from_correlation_function to construct an environment manually from one of these characteristic functions, or use a predefined sub-class such as the DrudeLorentzEnvironment, the UnderDampedEnvironment or the OhmicEnvironment.

Bosonic environments offer various ways to approximate the environment with a multi-exponential correlation function, which can be used for example in the HEOM solver. The approximated environment is represented as a ExponentialBosonicEnvironment.

All bosonic environments can be approximated by directly fitting their correlation function with a multi-exponential ansatz (approx_by_cf_fit) or by fitting their spectral density with a sum of Lorentzians (approx_by_sd_fit), which correspond to underdamped environments with known multi-exponential decompositions. Subclasses may offer additional approximation methods such as DrudeLorentzEnvironment.approx_by_matsubara or DrudeLorentzEnvironment.approx_by_pade in the case of a Drude-Lorentz environment.

Parameters:
Toptional, float

The temperature of this environment.

tagoptional, str, tuple or any other object

An identifier (name) for this environment.

approx_by_cf_fit(
tlist: ArrayLike,
target_rsme: float = 2e-05,
Nr_max: int = 10,
Ni_max: int = 10,
guess: list[float] = None,
lower: list[float] = None,
upper: list[float] = None,
sigma: float | ArrayLike = None,
maxfev: int = None,
full_ansatz: bool = False,
combine: bool = True,
tag: Any = None,
) tuple[ExponentialBosonicEnvironment, dict[str, Any]][source]

Generates an approximation to this environment by fitting its correlation function with a multi-exponential ansatz. The number of exponents is determined iteratively based on reducing the normalized root mean squared error below a given threshold.

Specifically, the real and imaginary parts are fit by the following model functions:

\[\begin{split}\operatorname{Re}[C(t)] = \sum_{k=1}^{N_r} \operatorname{Re}\Bigl[ (a_k + \mathrm i d_k) \mathrm e^{(b_k + \mathrm i c_k) t}\Bigl] , \\ \operatorname{Im}[C(t)] = \sum_{k=1}^{N_i} \operatorname{Im}\Bigl[ (a'_k + \mathrm i d'_k) \mathrm e^{(b'_k + \mathrm i c'_k) t} \Bigr].\end{split}\]

If the parameter full_ansatz is False, \(d_k\) and \(d'_k\) are set to zero and the model functions simplify to

\[\begin{split}\operatorname{Re}[C(t)] = \sum_{k=1}^{N_r} a_k e^{b_k t} \cos(c_{k} t) , \\ \operatorname{Im}[C(t)] = \sum_{k=1}^{N_i} a'_k e^{b'_k t} \sin(c'_{k} t) .\end{split}\]

The simplified version offers faster fits, however it fails for anomalous spectral densities with \(\operatorname{Im}[C(0)] \neq 0\) as \(\sin(0) = 0\).

Parameters:
tlistarray_like

The time range on which to perform the fit.

target_rmseoptional, float

Desired normalized root mean squared error (default 2e-5). Can be set to None to perform only one fit using the maximum number of modes (Nr_max, Ni_max).

Nr_maxoptional, int

The maximum number of modes to use for the fit of the real part (default 10).

Ni_maxoptional, int

The maximum number of modes to use for the fit of the imaginary part (default 10).

guessoptional, list of float

Initial guesses for the parameters \(a_k\), \(b_k\), etc. The same initial guesses are used for all values of k, and for the real and imaginary parts. If full_ansatz is True, guess is a list of size 4, otherwise, it is a list of size 3. If none of guess, lower and upper are provided, these parameters will be chosen automatically.

loweroptional, list of float

Lower bounds for the parameters \(a_k\), \(b_k\), etc. The same lower bounds are used for all values of k, and for the real and imaginary parts. If full_ansatz is True, lower is a list of size 4, otherwise, it is a list of size 3. If none of guess, lower and upper are provided, these parameters will be chosen automatically.

upperoptional, list of float

Upper bounds for the parameters \(a_k\), \(b_k\), etc. The same upper bounds are used for all values of k, and for the real and imaginary parts. If full_ansatz is True, upper is a list of size 4, otherwise, it is a list of size 3. If none of guess, lower and upper are provided, these parameters will be chosen automatically.

sigmaoptional, float or list of float

Adds an uncertainty to the correlation function of the environment, i.e., adds a leeway to the fit. This parameter is useful to adjust if the correlation function is very small in parts of the time range. For more details, see the documentation of scipy.optimize.curve_fit.

maxfevoptional, int

Number of times the parameters of the fit are allowed to vary during the optimization (per fit).

full_ansatzoptional, bool (default False)

If this is set to False, the parameters \(d_k\) are all set to zero. The full ansatz, including \(d_k\), usually leads to significantly slower fits, and some manual tuning of the guesses, lower and upper is usually needed. On the other hand, the full ansatz can lead to better fits with fewer exponents, especially for anomalous spectral densities with \(\operatorname{Im}[C(0)] \neq 0\) for which the simplified ansatz will always give \(\operatorname{Im}[C(0)] = 0\). When using the full ansatz with default values for the guesses and bounds, if the fit takes too long, we recommend choosing guesses and bounds manually.

combineoptional, bool (default True)

Whether to combine exponents with the same frequency. See combine for details.

tagoptional, str, tuple or any other object

An identifier (name) for the approximated environment. If not provided, a tag will be generated from the tag of this environment.

Returns:
approx_envExponentialBosonicEnvironment

The approximated environment with multi-exponential correlation function.

fit_infodictionary

A dictionary containing the following information about the fit.

“Nr”

The number of terms used to fit the real part of the correlation function.

“Ni”

The number of terms used to fit the imaginary part of the correlation function.

“fit_time_real”

The time the fit of the real part of the correlation function took in seconds.

“fit_time_imag”

The time the fit of the imaginary part of the correlation function took in seconds.

“rmse_real”

Normalized mean squared error obtained in the fit of the real part of the correlation function.

“rmse_imag”

Normalized mean squared error obtained in the fit of the imaginary part of the correlation function.

“params_real”

The fitted parameters (array of shape Nx3 or Nx4) for the real part of the correlation function.

“params_imag”

The fitted parameters (array of shape Nx3 or Nx4) for the imaginary part of the correlation function.

“summary”

A string that summarizes the information about the fit.

approx_by_sd_fit(
wlist: ArrayLike,
Nk: int = 1,
target_rmse: float = 5e-06,
Nmax: int = 10,
guess: list[float] = None,
lower: list[float] = None,
upper: list[float] = None,
sigma: float | ArrayLike = None,
maxfev: int = None,
combine: bool = True,
tag: Any = None,
) tuple[ExponentialBosonicEnvironment, dict[str, Any]][source]

Generates an approximation to this environment by fitting its spectral density with a sum of underdamped terms. Each underdamped term effectively acts like an underdamped environment. We use the known exponential decomposition of the underdamped environment, keeping Nk Matsubara terms for each. The number of underdamped terms is determined iteratively based on reducing the normalized root mean squared error below a given threshold.

Specifically, the spectral density is fit by the following model function:

\[J(\omega) = \sum_{k=1}^{N} \frac{2 a_k b_k \omega}{\left(\left( \omega + c_k \right)^2 + b_k^2 \right) \left(\left( \omega - c_k \right)^2 + b_k^2 \right)}\]
Parameters:
wlistarray_like

The frequency range on which to perform the fit.

Nkoptional, int

The number of Matsubara terms to keep in each mode (default 1).

target_rmseoptional, float

Desired normalized root mean squared error (default 5e-6). Can be set to None to perform only one fit using the maximum number of modes (Nmax).

Nmaxoptional, int

The maximum number of modes to use for the fit (default 10).

guessoptional, list of float

Initial guesses for the parameters \(a_k\), \(b_k\) and \(c_k\). The same initial guesses are used for all values of k. If none of guess, lower and upper are provided, these parameters will be chosen automatically.

loweroptional, list of float

Lower bounds for the parameters \(a_k\), \(b_k\) and \(c_k\). The same lower bounds are used for all values of k. If none of guess, lower and upper are provided, these parameters will be chosen automatically.

upperoptional, list of float

Upper bounds for the parameters \(a_k\), \(b_k\) and \(c_k\). The same upper bounds are used for all values of k. If none of guess, lower and upper are provided, these parameters will be chosen automatically.

sigmaoptional, float or list of float

Adds an uncertainty to the spectral density of the environment, i.e., adds a leeway to the fit. This parameter is useful to adjust if the spectral density is very small in parts of the frequency range. For more details, see the documentation of scipy.optimize.curve_fit.

maxfevoptional, int

Number of times the parameters of the fit are allowed to vary during the optimization (per fit).

combineoptional, bool (default True)

Whether to combine exponents with the same frequency. See combine for details.

tagoptional, str, tuple or any other object

An identifier (name) for the approximated environment. If not provided, a tag will be generated from the tag of this environment.

Returns:
approx_envExponentialBosonicEnvironment

The approximated environment with multi-exponential correlation function.

fit_infodictionary

A dictionary containing the following information about the fit.

“N”

The number of underdamped terms used in the fit.

“Nk”

The number of Matsubara modes included per underdamped term.

“fit_time”

The time the fit took in seconds.

“rmse”

Normalized mean squared error obtained in the fit.

“params”

The fitted parameters (array of shape Nx3).

“summary”

A string that summarizes the information about the fit.

abstract correlation_function(
t: float | ArrayLike,
*,
eps: float = 1e-10,
) float | ArrayLike[source]

The two-time auto-correlation function of this environment. See the Users Guide on bosonic environments for specifics on the definitions used by QuTiP.

If no analytical expression for the correlation function is known, it will be derived from the power spectrum via a fast fourier transform.

If no analytical expression for the power spectrum is known either, it will be derived from the spectral density. In this case, the temperature of this environment must be specified.

Parameters:
tarray_like or float

The times at which to evaluate the correlation function.

epsoptional, float

Used in case the power spectrum is derived from the spectral density; see the documentation of BosonicEnvironment.power_spectrum.

classmethod from_correlation_function(
C: Callable[[float], complex] | ArrayLike,
tlist: ArrayLike = None,
tMax: float = None,
*,
T: float = None,
tag: Any = None,
args: dict[str, Any] = None,
) BosonicEnvironment[source]

Constructs a bosonic environment from the provided correlation function. The provided function will only be used for times \(t \geq 0\). At times \(t < 0\), the symmetry relation \(C(-t) = C(t)^\ast\) is enforced.

Parameters:
Ccallable or array_like

The correlation function. Can be provided as a Python function or as an array. When using a function, the signature should be

C(t: array_like, **args) -> array_like

where t is time and args is a dict containing the other parameters of the function.

tlistoptional, array_like

The times where the correlation function is sampled (if it is provided as an array).

tMaxoptional, float

Specifies that the correlation function is essentially zero outside the interval [-tMax, tMax]. Used for numerical integration purposes.

Toptional, float

Environment temperature. (The spectral density of this environment can only be calculated from the correlation function if a temperature is provided.)

tagoptional, str, tuple or any other object

An identifier (name) for this environment.

argsoptional, dict

Extra arguments for the correlation function C.

classmethod from_power_spectrum(
S: Callable[[float], float] | ArrayLike,
wlist: ArrayLike = None,
wMax: float = None,
*,
T: float = None,
tag: Any = None,
args: dict[str, Any] = None,
) BosonicEnvironment[source]

Constructs a bosonic environment with the provided power spectrum.

Parameters:
Scallable or array_like

The power spectrum. Can be provided as a Python function or as an array. When using a function, the signature should be

S(w: array_like, **args) -> array_like

where w is the frequency and args is a dict containing the other parameters of the function.

wlistoptional, array_like

The frequencies where the power spectrum is sampled (if it is provided as an array).

wMaxoptional, float

Specifies that the power spectrum is essentially zero outside the interval [-wMax, wMax]. Used for numerical integration purposes.

Toptional, float

Environment temperature. (The spectral density of this environment can only be calculated from the powr spectrum if a temperature is provided.)

tagoptional, str, tuple or any other object

An identifier (name) for this environment.

argsoptional, dict

Extra arguments for the power spectrum S.

classmethod from_spectral_density(
J: Callable[[float], float] | ArrayLike,
wlist: ArrayLike = None,
wMax: float = None,
*,
T: float = None,
tag: Any = None,
args: dict[str, Any] = None,
) BosonicEnvironment[source]

Constructs a bosonic environment with the provided spectral density. The provided function will only be used for frequencies \(\omega > 0\). At frequencies \(\omega \leq 0\), the spectral density is zero according to the definition used by QuTiP. See the Users Guide on bosonic environments for a note on spectral densities with support at negative frequencies.

Parameters:
Jcallable or array_like

The spectral density. Can be provided as a Python function or as an array. When using a function, the signature should be

J(w: array_like, **args) -> array_like

where w is the frequency and args is a tuple containing the other parameters of the function.

wlistoptional, array_like

The frequencies where the spectral density is sampled (if it is provided as an array).

wMaxoptional, float

Specifies that the spectral density is essentially zero outside the interval [-wMax, wMax]. Used for numerical integration purposes.

Toptional, float

Environment temperature. (The correlation function and the power spectrum of this environment can only be calculated from the spectral density if a temperature is provided.)

tagoptional, str, tuple or any other object

An identifier (name) for this environment.

argsoptional, dict

Extra arguments for the spectral density S.

abstract power_spectrum(
w: float | ArrayLike,
*,
eps: float = 1e-10,
) float | ArrayLike[source]

The power spectrum of this environment. See the Users Guide on bosonic environments for specifics on the definitions used by QuTiP.

If no analytical expression for the power spectrum is known, it will be derived from the spectral density. In this case, the temperature of this environment must be specified.

If no analytical expression for the spectral density is known either, the power spectrum will instead be derived from the correlation function via a fast fourier transform.

Parameters:
warray_like or float

The frequencies at which to evaluate the power spectrum.

epsoptional, float

To derive the zero-frequency power spectrum from the spectral density, the spectral density must be differentiated numerically. In that case, this parameter is used as the finite difference in the numerical differentiation.

abstract spectral_density(w: float | ArrayLike) float | ArrayLike[source]

The spectral density of this environment. For negative frequencies, a value of zero will be returned. See the Users Guide on bosonic environments for specifics on the definitions used by QuTiP.

If no analytical expression for the spectral density is known, it will be derived from the power spectrum. In this case, the temperature of this environment must be specified.

If no analytical expression for the power spectrum is known either, it will be derived from the correlation function via a fast fourier transform.

Parameters:
warray_like or float

The frequencies at which to evaluate the spectral density.

class DrudeLorentzEnvironment(
T: float,
lam: float,
gamma: float,
*,
Nk: int = 10,
tag: Any = None,
)[source]

Bases: BosonicEnvironment

Describes a Drude-Lorentz bosonic environment with the spectral density

\[J(\omega) = \frac{2 \lambda \gamma \omega}{\gamma^{2}+\omega^{2}}\]

(see Eq. 15 in [BoFiN23]).

Parameters:
Tfloat

Environment temperature.

lamfloat

Coupling strength.

gammafloat

Spectral density cutoff frequency.

Nkoptional, int, defaults to 10

The number of Pade exponents to be used for the calculation of the correlation function.

tagoptional, str, tuple or any other object

An identifier (name) for this environment.

approx_by_cf_fit(
tlist: ArrayLike,
target_rsme: float = 2e-05,
Nr_max: int = 10,
Ni_max: int = 10,
guess: list[float] = None,
lower: list[float] = None,
upper: list[float] = None,
sigma: float | ArrayLike = None,
maxfev: int = None,
full_ansatz: bool = False,
combine: bool = True,
tag: Any = None,
) tuple[ExponentialBosonicEnvironment, dict[str, Any]]

Generates an approximation to this environment by fitting its correlation function with a multi-exponential ansatz. The number of exponents is determined iteratively based on reducing the normalized root mean squared error below a given threshold.

Specifically, the real and imaginary parts are fit by the following model functions:

\[\begin{split}\operatorname{Re}[C(t)] = \sum_{k=1}^{N_r} \operatorname{Re}\Bigl[ (a_k + \mathrm i d_k) \mathrm e^{(b_k + \mathrm i c_k) t}\Bigl] , \\ \operatorname{Im}[C(t)] = \sum_{k=1}^{N_i} \operatorname{Im}\Bigl[ (a'_k + \mathrm i d'_k) \mathrm e^{(b'_k + \mathrm i c'_k) t} \Bigr].\end{split}\]

If the parameter full_ansatz is False, \(d_k\) and \(d'_k\) are set to zero and the model functions simplify to

\[\begin{split}\operatorname{Re}[C(t)] = \sum_{k=1}^{N_r} a_k e^{b_k t} \cos(c_{k} t) , \\ \operatorname{Im}[C(t)] = \sum_{k=1}^{N_i} a'_k e^{b'_k t} \sin(c'_{k} t) .\end{split}\]

The simplified version offers faster fits, however it fails for anomalous spectral densities with \(\operatorname{Im}[C(0)] \neq 0\) as \(\sin(0) = 0\).

Parameters:
tlistarray_like

The time range on which to perform the fit.

target_rmseoptional, float

Desired normalized root mean squared error (default 2e-5). Can be set to None to perform only one fit using the maximum number of modes (Nr_max, Ni_max).

Nr_maxoptional, int

The maximum number of modes to use for the fit of the real part (default 10).

Ni_maxoptional, int

The maximum number of modes to use for the fit of the imaginary part (default 10).

guessoptional, list of float

Initial guesses for the parameters \(a_k\), \(b_k\), etc. The same initial guesses are used for all values of k, and for the real and imaginary parts. If full_ansatz is True, guess is a list of size 4, otherwise, it is a list of size 3. If none of guess, lower and upper are provided, these parameters will be chosen automatically.

loweroptional, list of float

Lower bounds for the parameters \(a_k\), \(b_k\), etc. The same lower bounds are used for all values of k, and for the real and imaginary parts. If full_ansatz is True, lower is a list of size 4, otherwise, it is a list of size 3. If none of guess, lower and upper are provided, these parameters will be chosen automatically.

upperoptional, list of float

Upper bounds for the parameters \(a_k\), \(b_k\), etc. The same upper bounds are used for all values of k, and for the real and imaginary parts. If full_ansatz is True, upper is a list of size 4, otherwise, it is a list of size 3. If none of guess, lower and upper are provided, these parameters will be chosen automatically.

sigmaoptional, float or list of float

Adds an uncertainty to the correlation function of the environment, i.e., adds a leeway to the fit. This parameter is useful to adjust if the correlation function is very small in parts of the time range. For more details, see the documentation of scipy.optimize.curve_fit.

maxfevoptional, int

Number of times the parameters of the fit are allowed to vary during the optimization (per fit).

full_ansatzoptional, bool (default False)

If this is set to False, the parameters \(d_k\) are all set to zero. The full ansatz, including \(d_k\), usually leads to significantly slower fits, and some manual tuning of the guesses, lower and upper is usually needed. On the other hand, the full ansatz can lead to better fits with fewer exponents, especially for anomalous spectral densities with \(\operatorname{Im}[C(0)] \neq 0\) for which the simplified ansatz will always give \(\operatorname{Im}[C(0)] = 0\). When using the full ansatz with default values for the guesses and bounds, if the fit takes too long, we recommend choosing guesses and bounds manually.

combineoptional, bool (default True)

Whether to combine exponents with the same frequency. See combine for details.

tagoptional, str, tuple or any other object

An identifier (name) for the approximated environment. If not provided, a tag will be generated from the tag of this environment.

Returns:
approx_envExponentialBosonicEnvironment

The approximated environment with multi-exponential correlation function.

fit_infodictionary

A dictionary containing the following information about the fit.

“Nr”

The number of terms used to fit the real part of the correlation function.

“Ni”

The number of terms used to fit the imaginary part of the correlation function.

“fit_time_real”

The time the fit of the real part of the correlation function took in seconds.

“fit_time_imag”

The time the fit of the imaginary part of the correlation function took in seconds.

“rmse_real”

Normalized mean squared error obtained in the fit of the real part of the correlation function.

“rmse_imag”

Normalized mean squared error obtained in the fit of the imaginary part of the correlation function.

“params_real”

The fitted parameters (array of shape Nx3 or Nx4) for the real part of the correlation function.

“params_imag”

The fitted parameters (array of shape Nx3 or Nx4) for the imaginary part of the correlation function.

“summary”

A string that summarizes the information about the fit.

approx_by_matsubara(
Nk: int,
combine: bool = True,
compute_delta: Literal[False] = False,
tag: Any = None,
) ExponentialBosonicEnvironment[source]
approx_by_matsubara(
Nk: int,
combine: bool = True,
compute_delta: Literal[True] = True,
tag: Any = None,
) tuple[ExponentialBosonicEnvironment, float]

Generates an approximation to this environment by truncating its Matsubara expansion.

Parameters:
Nkint

Number of Matsubara terms to include. In total, the real part of the correlation function will include Nk+1 terms and the imaginary part 1 term.

combinebool, default True

Whether to combine exponents with the same frequency.

compute_deltabool, default False

Whether to compute and return the approximation discrepancy (see below).

tagoptional, str, tuple or any other object

An identifier (name) for the approximated environment. If not provided, a tag will be generated from the tag of this environment.

Returns:
approx_envExponentialBosonicEnvironment

The approximated environment with multi-exponential correlation function.

deltafloat, optional

The approximation discrepancy. That is, the difference between the true correlation function of the Drude-Lorentz environment and the sum of the Nk exponential terms is approximately 2 * delta * dirac(t), where dirac(t) denotes the Dirac delta function. It can be used to create a “terminator” term to add to the system dynamics to take this discrepancy into account, see system_terminator.

approx_by_pade(
Nk: int,
combine: bool = True,
compute_delta: Literal[False] = False,
tag: Any = None,
) ExponentialBosonicEnvironment[source]
approx_by_pade(
Nk: int,
combine: bool = True,
compute_delta: Literal[True] = True,
tag: Any = None,
) tuple[ExponentialBosonicEnvironment, float]

Generates an approximation to this environment by truncating its Pade expansion.

Parameters:
Nkint

Number of Pade terms to include. In total, the real part of the correlation function will include Nk+1 terms and the imaginary part 1 term.

combinebool, default True

Whether to combine exponents with the same frequency.

compute_deltabool, default False

Whether to compute and return the approximation discrepancy (see below).

tagoptional, str, tuple or any other object

An identifier (name) for the approximated environment. If not provided, a tag will be generated from the tag of this environment.

Returns:
approx_envExponentialBosonicEnvironment

The approximated environment with multi-exponential correlation function.

deltafloat, optional

The approximation discrepancy. That is, the difference between the true correlation function of the Drude-Lorentz environment and the sum of the Nk exponential terms is approximately 2 * delta * dirac(t), where dirac(t) denotes the Dirac delta function. It can be used to create a “terminator” term to add to the system dynamics to take this discrepancy into account, see system_terminator.

approx_by_sd_fit(
wlist: ArrayLike,
Nk: int = 1,
target_rmse: float = 5e-06,
Nmax: int = 10,
guess: list[float] = None,
lower: list[float] = None,
upper: list[float] = None,
sigma: float | ArrayLike = None,
maxfev: int = None,
combine: bool = True,
tag: Any = None,
) tuple[ExponentialBosonicEnvironment, dict[str, Any]]

Generates an approximation to this environment by fitting its spectral density with a sum of underdamped terms. Each underdamped term effectively acts like an underdamped environment. We use the known exponential decomposition of the underdamped environment, keeping Nk Matsubara terms for each. The number of underdamped terms is determined iteratively based on reducing the normalized root mean squared error below a given threshold.

Specifically, the spectral density is fit by the following model function:

\[J(\omega) = \sum_{k=1}^{N} \frac{2 a_k b_k \omega}{\left(\left( \omega + c_k \right)^2 + b_k^2 \right) \left(\left( \omega - c_k \right)^2 + b_k^2 \right)}\]
Parameters:
wlistarray_like

The frequency range on which to perform the fit.

Nkoptional, int

The number of Matsubara terms to keep in each mode (default 1).

target_rmseoptional, float

Desired normalized root mean squared error (default 5e-6). Can be set to None to perform only one fit using the maximum number of modes (Nmax).

Nmaxoptional, int

The maximum number of modes to use for the fit (default 10).

guessoptional, list of float

Initial guesses for the parameters \(a_k\), \(b_k\) and \(c_k\). The same initial guesses are used for all values of k. If none of guess, lower and upper are provided, these parameters will be chosen automatically.

loweroptional, list of float

Lower bounds for the parameters \(a_k\), \(b_k\) and \(c_k\). The same lower bounds are used for all values of k. If none of guess, lower and upper are provided, these parameters will be chosen automatically.

upperoptional, list of float

Upper bounds for the parameters \(a_k\), \(b_k\) and \(c_k\). The same upper bounds are used for all values of k. If none of guess, lower and upper are provided, these parameters will be chosen automatically.

sigmaoptional, float or list of float

Adds an uncertainty to the spectral density of the environment, i.e., adds a leeway to the fit. This parameter is useful to adjust if the spectral density is very small in parts of the frequency range. For more details, see the documentation of scipy.optimize.curve_fit.

maxfevoptional, int

Number of times the parameters of the fit are allowed to vary during the optimization (per fit).

combineoptional, bool (default True)

Whether to combine exponents with the same frequency. See combine for details.

tagoptional, str, tuple or any other object

An identifier (name) for the approximated environment. If not provided, a tag will be generated from the tag of this environment.

Returns:
approx_envExponentialBosonicEnvironment

The approximated environment with multi-exponential correlation function.

fit_infodictionary

A dictionary containing the following information about the fit.

“N”

The number of underdamped terms used in the fit.

“Nk”

The number of Matsubara modes included per underdamped term.

“fit_time”

The time the fit took in seconds.

“rmse”

Normalized mean squared error obtained in the fit.

“params”

The fitted parameters (array of shape Nx3).

“summary”

A string that summarizes the information about the fit.

correlation_function(
t: float | ArrayLike,
Nk: int = None,
**kwargs,
) float | ArrayLike[source]

Calculates the two-time auto-correlation function of the Drude-Lorentz environment. The calculation is performed by summing a large number of exponents of the Pade expansion.

Parameters:
tarray_like or float

The time at which to evaluate the correlation function.

Nkint, optional

The number of exponents to use. If not provided, then the value that was provided when the class was instantiated is used.

power_spectrum(
w: float | ArrayLike,
**kwargs,
) float | ArrayLike[source]

Calculates the power spectrum of the Drude-Lorentz environment.

Parameters:
warray_like or float

The frequency at which to evaluate the power spectrum.

spectral_density(
w: float | ArrayLike,
) float | ArrayLike[source]

Calculates the Drude-Lorentz spectral density.

Parameters:
warray_like or float

Energy of the mode.

class UnderDampedEnvironment(
T: float,
lam: float,
gamma: float,
w0: float,
*,
tag: Any = None,
)[source]

Bases: BosonicEnvironment

Describes an underdamped environment with the spectral density

\[J(\omega) = \frac{\lambda^{2} \Gamma \omega}{(\omega_0^{2}- \omega^{2})^{2}+ \Gamma^{2} \omega^{2}}\]

(see Eq. 16 in [BoFiN23]).

Parameters:
Tfloat

Environment temperature.

lamfloat

Coupling strength.

gammafloat

Spectral density cutoff frequency.

w0float

Spectral density resonance frequency.

tagoptional, str, tuple or any other object

An identifier (name) for this environment.

approx_by_cf_fit(
tlist: ArrayLike,
target_rsme: float = 2e-05,
Nr_max: int = 10,
Ni_max: int = 10,
guess: list[float] = None,
lower: list[float] = None,
upper: list[float] = None,
sigma: float | ArrayLike = None,
maxfev: int = None,
full_ansatz: bool = False,
combine: bool = True,
tag: Any = None,
) tuple[ExponentialBosonicEnvironment, dict[str, Any]]

Generates an approximation to this environment by fitting its correlation function with a multi-exponential ansatz. The number of exponents is determined iteratively based on reducing the normalized root mean squared error below a given threshold.

Specifically, the real and imaginary parts are fit by the following model functions:

\[\begin{split}\operatorname{Re}[C(t)] = \sum_{k=1}^{N_r} \operatorname{Re}\Bigl[ (a_k + \mathrm i d_k) \mathrm e^{(b_k + \mathrm i c_k) t}\Bigl] , \\ \operatorname{Im}[C(t)] = \sum_{k=1}^{N_i} \operatorname{Im}\Bigl[ (a'_k + \mathrm i d'_k) \mathrm e^{(b'_k + \mathrm i c'_k) t} \Bigr].\end{split}\]

If the parameter full_ansatz is False, \(d_k\) and \(d'_k\) are set to zero and the model functions simplify to

\[\begin{split}\operatorname{Re}[C(t)] = \sum_{k=1}^{N_r} a_k e^{b_k t} \cos(c_{k} t) , \\ \operatorname{Im}[C(t)] = \sum_{k=1}^{N_i} a'_k e^{b'_k t} \sin(c'_{k} t) .\end{split}\]

The simplified version offers faster fits, however it fails for anomalous spectral densities with \(\operatorname{Im}[C(0)] \neq 0\) as \(\sin(0) = 0\).

Parameters:
tlistarray_like

The time range on which to perform the fit.

target_rmseoptional, float

Desired normalized root mean squared error (default 2e-5). Can be set to None to perform only one fit using the maximum number of modes (Nr_max, Ni_max).

Nr_maxoptional, int

The maximum number of modes to use for the fit of the real part (default 10).

Ni_maxoptional, int

The maximum number of modes to use for the fit of the imaginary part (default 10).

guessoptional, list of float

Initial guesses for the parameters \(a_k\), \(b_k\), etc. The same initial guesses are used for all values of k, and for the real and imaginary parts. If full_ansatz is True, guess is a list of size 4, otherwise, it is a list of size 3. If none of guess, lower and upper are provided, these parameters will be chosen automatically.

loweroptional, list of float

Lower bounds for the parameters \(a_k\), \(b_k\), etc. The same lower bounds are used for all values of k, and for the real and imaginary parts. If full_ansatz is True, lower is a list of size 4, otherwise, it is a list of size 3. If none of guess, lower and upper are provided, these parameters will be chosen automatically.

upperoptional, list of float

Upper bounds for the parameters \(a_k\), \(b_k\), etc. The same upper bounds are used for all values of k, and for the real and imaginary parts. If full_ansatz is True, upper is a list of size 4, otherwise, it is a list of size 3. If none of guess, lower and upper are provided, these parameters will be chosen automatically.

sigmaoptional, float or list of float

Adds an uncertainty to the correlation function of the environment, i.e., adds a leeway to the fit. This parameter is useful to adjust if the correlation function is very small in parts of the time range. For more details, see the documentation of scipy.optimize.curve_fit.

maxfevoptional, int

Number of times the parameters of the fit are allowed to vary during the optimization (per fit).

full_ansatzoptional, bool (default False)

If this is set to False, the parameters \(d_k\) are all set to zero. The full ansatz, including \(d_k\), usually leads to significantly slower fits, and some manual tuning of the guesses, lower and upper is usually needed. On the other hand, the full ansatz can lead to better fits with fewer exponents, especially for anomalous spectral densities with \(\operatorname{Im}[C(0)] \neq 0\) for which the simplified ansatz will always give \(\operatorname{Im}[C(0)] = 0\). When using the full ansatz with default values for the guesses and bounds, if the fit takes too long, we recommend choosing guesses and bounds manually.

combineoptional, bool (default True)

Whether to combine exponents with the same frequency. See combine for details.

tagoptional, str, tuple or any other object

An identifier (name) for the approximated environment. If not provided, a tag will be generated from the tag of this environment.

Returns:
approx_envExponentialBosonicEnvironment

The approximated environment with multi-exponential correlation function.

fit_infodictionary

A dictionary containing the following information about the fit.

“Nr”

The number of terms used to fit the real part of the correlation function.

“Ni”

The number of terms used to fit the imaginary part of the correlation function.

“fit_time_real”

The time the fit of the real part of the correlation function took in seconds.

“fit_time_imag”

The time the fit of the imaginary part of the correlation function took in seconds.

“rmse_real”

Normalized mean squared error obtained in the fit of the real part of the correlation function.

“rmse_imag”

Normalized mean squared error obtained in the fit of the imaginary part of the correlation function.

“params_real”

The fitted parameters (array of shape Nx3 or Nx4) for the real part of the correlation function.

“params_imag”

The fitted parameters (array of shape Nx3 or Nx4) for the imaginary part of the correlation function.

“summary”

A string that summarizes the information about the fit.

approx_by_matsubara(
Nk: int,
combine: bool = True,
tag: Any = None,
) ExponentialBosonicEnvironment[source]

Generates an approximation to this environment by truncating its Matsubara expansion.

Parameters:
Nkint

Number of Matsubara terms to include. In total, the real part of the correlation function will include Nk+2 terms and the imaginary part 2 terms.

combinebool, default True

Whether to combine exponents with the same frequency.

tagoptional, str, tuple or any other object

An identifier (name) for the approximated environment. If not provided, a tag will be generated from the tag of this environment.

Returns:
ExponentialBosonicEnvironment

The approximated environment with multi-exponential correlation function.

approx_by_sd_fit(
wlist: ArrayLike,
Nk: int = 1,
target_rmse: float = 5e-06,
Nmax: int = 10,
guess: list[float] = None,
lower: list[float] = None,
upper: list[float] = None,
sigma: float | ArrayLike = None,
maxfev: int = None,
combine: bool = True,
tag: Any = None,
) tuple[ExponentialBosonicEnvironment, dict[str, Any]]

Generates an approximation to this environment by fitting its spectral density with a sum of underdamped terms. Each underdamped term effectively acts like an underdamped environment. We use the known exponential decomposition of the underdamped environment, keeping Nk Matsubara terms for each. The number of underdamped terms is determined iteratively based on reducing the normalized root mean squared error below a given threshold.

Specifically, the spectral density is fit by the following model function:

\[J(\omega) = \sum_{k=1}^{N} \frac{2 a_k b_k \omega}{\left(\left( \omega + c_k \right)^2 + b_k^2 \right) \left(\left( \omega - c_k \right)^2 + b_k^2 \right)}\]
Parameters:
wlistarray_like

The frequency range on which to perform the fit.

Nkoptional, int

The number of Matsubara terms to keep in each mode (default 1).

target_rmseoptional, float

Desired normalized root mean squared error (default 5e-6). Can be set to None to perform only one fit using the maximum number of modes (Nmax).

Nmaxoptional, int

The maximum number of modes to use for the fit (default 10).

guessoptional, list of float

Initial guesses for the parameters \(a_k\), \(b_k\) and \(c_k\). The same initial guesses are used for all values of k. If none of guess, lower and upper are provided, these parameters will be chosen automatically.

loweroptional, list of float

Lower bounds for the parameters \(a_k\), \(b_k\) and \(c_k\). The same lower bounds are used for all values of k. If none of guess, lower and upper are provided, these parameters will be chosen automatically.

upperoptional, list of float

Upper bounds for the parameters \(a_k\), \(b_k\) and \(c_k\). The same upper bounds are used for all values of k. If none of guess, lower and upper are provided, these parameters will be chosen automatically.

sigmaoptional, float or list of float

Adds an uncertainty to the spectral density of the environment, i.e., adds a leeway to the fit. This parameter is useful to adjust if the spectral density is very small in parts of the frequency range. For more details, see the documentation of scipy.optimize.curve_fit.

maxfevoptional, int

Number of times the parameters of the fit are allowed to vary during the optimization (per fit).

combineoptional, bool (default True)

Whether to combine exponents with the same frequency. See combine for details.

tagoptional, str, tuple or any other object

An identifier (name) for the approximated environment. If not provided, a tag will be generated from the tag of this environment.

Returns:
approx_envExponentialBosonicEnvironment

The approximated environment with multi-exponential correlation function.

fit_infodictionary

A dictionary containing the following information about the fit.

“N”

The number of underdamped terms used in the fit.

“Nk”

The number of Matsubara modes included per underdamped term.

“fit_time”

The time the fit took in seconds.

“rmse”

Normalized mean squared error obtained in the fit.

“params”

The fitted parameters (array of shape Nx3).

“summary”

A string that summarizes the information about the fit.

correlation_function(
t: float | ArrayLike,
**kwargs,
) float | ArrayLike[source]

Calculates the two-time auto-correlation function of the underdamped environment.

Parameters:
tarray_like or float

The time at which to evaluate the correlation function.

power_spectrum(
w: float | ArrayLike,
**kwargs,
) float | ArrayLike[source]

Calculates the power spectrum of the underdamped environment.

Parameters:
warray_like or float

The frequency at which to evaluate the power spectrum.

spectral_density(
w: float | ArrayLike,
) float | ArrayLike[source]

Calculates the underdamped spectral density.

Parameters:
warray_like or float

Energy of the mode.

class OhmicEnvironment(
T: float,
alpha: float,
wc: float,
s: float,
*,
tag: Any = None,
)[source]

Bases: BosonicEnvironment

Describes Ohmic environments as well as sub- or super-Ohmic environments (depending on the choice of the parameter s). The spectral density is

\[J(\omega) = \alpha \frac{\omega^s}{\omega_c^{s-1}} e^{-\omega / \omega_c} .\]

This class requires the mpmath module to be installed.

Parameters:
Tfloat

Temperature of the environment.

alphafloat

Coupling strength.

wcfloat

Cutoff parameter.

sfloat

Power of omega in the spectral density.

tagoptional, str, tuple or any other object

An identifier (name) for this environment.

approx_by_cf_fit(
tlist: ArrayLike,
target_rsme: float = 2e-05,
Nr_max: int = 10,
Ni_max: int = 10,
guess: list[float] = None,
lower: list[float] = None,
upper: list[float] = None,
sigma: float | ArrayLike = None,
maxfev: int = None,
full_ansatz: bool = False,
combine: bool = True,
tag: Any = None,
) tuple[ExponentialBosonicEnvironment, dict[str, Any]]

Generates an approximation to this environment by fitting its correlation function with a multi-exponential ansatz. The number of exponents is determined iteratively based on reducing the normalized root mean squared error below a given threshold.

Specifically, the real and imaginary parts are fit by the following model functions:

\[\begin{split}\operatorname{Re}[C(t)] = \sum_{k=1}^{N_r} \operatorname{Re}\Bigl[ (a_k + \mathrm i d_k) \mathrm e^{(b_k + \mathrm i c_k) t}\Bigl] , \\ \operatorname{Im}[C(t)] = \sum_{k=1}^{N_i} \operatorname{Im}\Bigl[ (a'_k + \mathrm i d'_k) \mathrm e^{(b'_k + \mathrm i c'_k) t} \Bigr].\end{split}\]

If the parameter full_ansatz is False, \(d_k\) and \(d'_k\) are set to zero and the model functions simplify to

\[\begin{split}\operatorname{Re}[C(t)] = \sum_{k=1}^{N_r} a_k e^{b_k t} \cos(c_{k} t) , \\ \operatorname{Im}[C(t)] = \sum_{k=1}^{N_i} a'_k e^{b'_k t} \sin(c'_{k} t) .\end{split}\]

The simplified version offers faster fits, however it fails for anomalous spectral densities with \(\operatorname{Im}[C(0)] \neq 0\) as \(\sin(0) = 0\).

Parameters:
tlistarray_like

The time range on which to perform the fit.

target_rmseoptional, float

Desired normalized root mean squared error (default 2e-5). Can be set to None to perform only one fit using the maximum number of modes (Nr_max, Ni_max).

Nr_maxoptional, int

The maximum number of modes to use for the fit of the real part (default 10).

Ni_maxoptional, int

The maximum number of modes to use for the fit of the imaginary part (default 10).

guessoptional, list of float

Initial guesses for the parameters \(a_k\), \(b_k\), etc. The same initial guesses are used for all values of k, and for the real and imaginary parts. If full_ansatz is True, guess is a list of size 4, otherwise, it is a list of size 3. If none of guess, lower and upper are provided, these parameters will be chosen automatically.

loweroptional, list of float

Lower bounds for the parameters \(a_k\), \(b_k\), etc. The same lower bounds are used for all values of k, and for the real and imaginary parts. If full_ansatz is True, lower is a list of size 4, otherwise, it is a list of size 3. If none of guess, lower and upper are provided, these parameters will be chosen automatically.

upperoptional, list of float

Upper bounds for the parameters \(a_k\), \(b_k\), etc. The same upper bounds are used for all values of k, and for the real and imaginary parts. If full_ansatz is True, upper is a list of size 4, otherwise, it is a list of size 3. If none of guess, lower and upper are provided, these parameters will be chosen automatically.

sigmaoptional, float or list of float

Adds an uncertainty to the correlation function of the environment, i.e., adds a leeway to the fit. This parameter is useful to adjust if the correlation function is very small in parts of the time range. For more details, see the documentation of scipy.optimize.curve_fit.

maxfevoptional, int

Number of times the parameters of the fit are allowed to vary during the optimization (per fit).

full_ansatzoptional, bool (default False)

If this is set to False, the parameters \(d_k\) are all set to zero. The full ansatz, including \(d_k\), usually leads to significantly slower fits, and some manual tuning of the guesses, lower and upper is usually needed. On the other hand, the full ansatz can lead to better fits with fewer exponents, especially for anomalous spectral densities with \(\operatorname{Im}[C(0)] \neq 0\) for which the simplified ansatz will always give \(\operatorname{Im}[C(0)] = 0\). When using the full ansatz with default values for the guesses and bounds, if the fit takes too long, we recommend choosing guesses and bounds manually.

combineoptional, bool (default True)

Whether to combine exponents with the same frequency. See combine for details.

tagoptional, str, tuple or any other object

An identifier (name) for the approximated environment. If not provided, a tag will be generated from the tag of this environment.

Returns:
approx_envExponentialBosonicEnvironment

The approximated environment with multi-exponential correlation function.

fit_infodictionary

A dictionary containing the following information about the fit.

“Nr”

The number of terms used to fit the real part of the correlation function.

“Ni”

The number of terms used to fit the imaginary part of the correlation function.

“fit_time_real”

The time the fit of the real part of the correlation function took in seconds.

“fit_time_imag”

The time the fit of the imaginary part of the correlation function took in seconds.

“rmse_real”

Normalized mean squared error obtained in the fit of the real part of the correlation function.

“rmse_imag”

Normalized mean squared error obtained in the fit of the imaginary part of the correlation function.

“params_real”

The fitted parameters (array of shape Nx3 or Nx4) for the real part of the correlation function.

“params_imag”

The fitted parameters (array of shape Nx3 or Nx4) for the imaginary part of the correlation function.

“summary”

A string that summarizes the information about the fit.

approx_by_sd_fit(
wlist: ArrayLike,
Nk: int = 1,
target_rmse: float = 5e-06,
Nmax: int = 10,
guess: list[float] = None,
lower: list[float] = None,
upper: list[float] = None,
sigma: float | ArrayLike = None,
maxfev: int = None,
combine: bool = True,
tag: Any = None,
) tuple[ExponentialBosonicEnvironment, dict[str, Any]]

Generates an approximation to this environment by fitting its spectral density with a sum of underdamped terms. Each underdamped term effectively acts like an underdamped environment. We use the known exponential decomposition of the underdamped environment, keeping Nk Matsubara terms for each. The number of underdamped terms is determined iteratively based on reducing the normalized root mean squared error below a given threshold.

Specifically, the spectral density is fit by the following model function:

\[J(\omega) = \sum_{k=1}^{N} \frac{2 a_k b_k \omega}{\left(\left( \omega + c_k \right)^2 + b_k^2 \right) \left(\left( \omega - c_k \right)^2 + b_k^2 \right)}\]
Parameters:
wlistarray_like

The frequency range on which to perform the fit.

Nkoptional, int

The number of Matsubara terms to keep in each mode (default 1).

target_rmseoptional, float

Desired normalized root mean squared error (default 5e-6). Can be set to None to perform only one fit using the maximum number of modes (Nmax).

Nmaxoptional, int

The maximum number of modes to use for the fit (default 10).

guessoptional, list of float

Initial guesses for the parameters \(a_k\), \(b_k\) and \(c_k\). The same initial guesses are used for all values of k. If none of guess, lower and upper are provided, these parameters will be chosen automatically.

loweroptional, list of float

Lower bounds for the parameters \(a_k\), \(b_k\) and \(c_k\). The same lower bounds are used for all values of k. If none of guess, lower and upper are provided, these parameters will be chosen automatically.

upperoptional, list of float

Upper bounds for the parameters \(a_k\), \(b_k\) and \(c_k\). The same upper bounds are used for all values of k. If none of guess, lower and upper are provided, these parameters will be chosen automatically.

sigmaoptional, float or list of float

Adds an uncertainty to the spectral density of the environment, i.e., adds a leeway to the fit. This parameter is useful to adjust if the spectral density is very small in parts of the frequency range. For more details, see the documentation of scipy.optimize.curve_fit.

maxfevoptional, int

Number of times the parameters of the fit are allowed to vary during the optimization (per fit).

combineoptional, bool (default True)

Whether to combine exponents with the same frequency. See combine for details.

tagoptional, str, tuple or any other object

An identifier (name) for the approximated environment. If not provided, a tag will be generated from the tag of this environment.

Returns:
approx_envExponentialBosonicEnvironment

The approximated environment with multi-exponential correlation function.

fit_infodictionary

A dictionary containing the following information about the fit.

“N”

The number of underdamped terms used in the fit.

“Nk”

The number of Matsubara modes included per underdamped term.

“fit_time”

The time the fit took in seconds.

“rmse”

Normalized mean squared error obtained in the fit.

“params”

The fitted parameters (array of shape Nx3).

“summary”

A string that summarizes the information about the fit.

correlation_function(
t: float | ArrayLike,
**kwargs,
) float | ArrayLike[source]

Calculates the correlation function of an Ohmic environment using the formula

\[C(t)= \frac{1}{\pi} \alpha w_{c}^{1-s} \beta^{-(s+1)} \Gamma(s+1) \left[ \zeta\left(s+1,\frac{1+\beta w_{c} -i w_{c} t}{\beta w_{c}} \right) +\zeta\left(s+1,\frac{1+ i w_{c} t}{\beta w_{c}}\right) \right] ,\]

where \(\Gamma\) is the gamma function, and \(\zeta\) the Riemann zeta function.

Parameters:
tarray_like or float

The time at which to evaluate the correlation function.

power_spectrum(
w: float | ArrayLike,
**kwargs,
) float | ArrayLike[source]

Calculates the power spectrum of the Ohmic environment.

Parameters:
warray_like or float

The frequency at which to evaluate the power spectrum.

spectral_density(w: float | ArrayLike) float | ArrayLike[source]

Calculates the spectral density of the Ohmic environment.

Parameters:
warray_like or float

Energy of the mode.

class CFExponent(
type: str | CFExponent.ExponentType,
ck: complex,
vk: complex,
ck2: complex = None,
tag: Any = None,
)[source]

Represents a single exponent (naively, an excitation mode) within an exponential decomposition of the correlation function of a environment.

Parameters:
type{“R”, “I”, “RI”, “+”, “-”} or one of CFExponent.types

The type of exponent.

“R” and “I” are bosonic exponents that appear in the real and imaginary parts of the correlation expansion, respectively.

“RI” is a combined bosonic exponent that appears in both the real and imaginary parts of the correlation expansion. The combined exponent has a single vk. The ck is the coefficient in the real expansion and ck2 is the coefficient in the imaginary expansion.

“+” and “-” are fermionic exponents.

ckcomplex

The coefficient of the excitation term.

vkcomplex

The frequency of the exponent of the excitation term.

ck2optional, complex

For exponents of type “RI” this is the coefficient of the term in the imaginary expansion (and ck is the coefficient in the real expansion).

tagoptional, str, tuple or any other object

A label for the exponent (often the name of the environment). It defaults to None.

Attributes:
fermionicbool

True if the type of the exponent is a Fermionic type (i.e. either “+” or “-”) and False otherwise.

coefficientcomplex

The coefficient of this excitation term in the total correlation function (including real and imaginary part).

exponentcomplex

The frequency of the exponent of the excitation term. (Alias for vk.)

All of the parameters are also available as attributes.
types

alias of ExponentType

class ExponentialBosonicEnvironment(
ck_real: ArrayLike = None,
vk_real: ArrayLike = None,
ck_imag: ArrayLike = None,
vk_imag: ArrayLike = None,
*,
exponents: Sequence[CFExponent] = None,
combine: bool = True,
T: float = None,
tag: Any = None,
)[source]

Bases: BosonicEnvironment

Bosonic environment that is specified through an exponential decomposition of its correlation function. The list of coefficients and exponents in the decomposition may either be passed through the four lists ck_real, vk_real, ck_imag, vk_imag, or as a list of bosonic CFExponent objects.

Parameters:
ck_reallist of complex

The coefficients of the expansion terms for the real part of the correlation function. The corresponding frequencies are passed as vk_real.

vk_reallist of complex

The frequencies (exponents) of the expansion terms for the real part of the correlation function. The corresponding coefficients are passed as ck_real.

ck_imaglist of complex

The coefficients of the expansion terms in the imaginary part of the correlation function. The corresponding frequencies are passed as vk_imag.

vk_imaglist of complex

The frequencies (exponents) of the expansion terms for the imaginary part of the correlation function. The corresponding coefficients are passed as ck_imag.

exponentslist of CFExponent

The expansion coefficients and exponents of both the real and the imaginary parts of the correlation function as CFExponent objects.

combinebool, default True

Whether to combine exponents with the same frequency. See combine for details.

T: optional, float

The temperature of the environment.

tagoptional, str, tuple or any other object

An identifier (name) for this environment.

classmethod combine(
exponents: Sequence[CFExponent],
rtol: float = 1e-05,
atol: float = 1e-07,
) Sequence[CFExponent][source]

Group bosonic exponents with the same frequency and return a single exponent for each frequency present.

Parameters:
exponentslist of CFExponent

The list of exponents to combine.

rtolfloat, default 1e-5

The relative tolerance to use to when comparing frequencies.

atolfloat, default 1e-7

The absolute tolerance to use to when comparing frequencies.

Returns:
list of CFExponent

The new reduced list of exponents.

correlation_function(
t: float | ArrayLike,
**kwargs,
) float | ArrayLike[source]

Computes the correlation function represented by this exponential decomposition.

Parameters:
tarray_like or float

The time at which to evaluate the correlation function.

power_spectrum(
w: float | ArrayLike,
**kwargs,
) float | ArrayLike[source]

Calculates the power spectrum corresponding to the multi-exponential correlation function.

Parameters:
warray_like or float

The frequency at which to evaluate the power spectrum.

spectral_density(
w: float | ArrayLike,
) float | ArrayLike[source]

Calculates the spectral density corresponding to the multi-exponential correlation function.

Parameters:
warray_like or float

Energy of the mode.

system_terminator(
Q: Qobj,
delta: float,
) Qobj[source]

Constructs the terminator for a given approximation discrepancy.

Parameters:
QQobj

The system coupling operator.

deltafloat

The approximation discrepancy of approximating an environment with a finite number of exponentials, see for example DrudeLorentzEnvironment.approx_by_matsubara.

Returns:
terminatorQobj

A superoperator acting on the system Hilbert space. Liouvillian term representing the contribution to the system-environment dynamics of all neglected expansion terms. It should be used by adding it to the system Liouvillian (i.e. liouvillian(H_sys)).

Fermionic Environments

class FermionicEnvironment(T: float = None, mu: float = None, tag: Any = None)[source]

The fermionic environment of an open quantum system. It is characterized by its spectral density, temperature and chemical potential or, equivalently, by its power spectra or its two-time auto-correlation functions.

This class is included as a counterpart to BosonicEnvironment, but it currently does not support all features that the bosonic environment does. In particular, fermionic environments cannot be constructed from manually specified spectral densities, power spectra or correlation functions. The only types of fermionic environment implemented at this time are Lorentzian environments (LorentzianEnvironment) and environments with multi-exponential correlation functions (ExponentialFermionicEnvironment).

Parameters:
Toptional, float

The temperature of this environment.

muoptional, float

The chemical potential of this environment.

tagoptional, str, tuple or any other object

An identifier (name) for this environment.

abstract correlation_function_minus(
t: float | ArrayLike,
) float | ArrayLike[source]

The “-“-branch of the auto-correlation function of this environment. See the Users Guide on fermionic environments for specifics on the definitions used by QuTiP.

Parameters:
tarray_like or float

The times at which to evaluate the correlation function.

abstract correlation_function_plus(
t: float | ArrayLike,
) float | ArrayLike[source]

The “+”-branch of the auto-correlation function of this environment. See the Users Guide on fermionic environments for specifics on the definitions used by QuTiP.

Parameters:
tarray_like or float

The times at which to evaluate the correlation function.

classmethod from_correlation_functions(
**kwargs,
) FermionicEnvironment[source]

User-defined fermionic environments are currently not implemented.

classmethod from_power_spectra(
**kwargs,
) FermionicEnvironment[source]

User-defined fermionic environments are currently not implemented.

abstract power_spectrum_minus(
w: float | ArrayLike,
) float | ArrayLike[source]

The “-“-branch of the power spectrum of this environment. See the Users Guide on fermionic environments for specifics on the definitions used by QuTiP.

Parameters:
warray_like or float

The frequencies at which to evaluate the power spectrum.

abstract power_spectrum_plus(
w: float | ArrayLike,
) float | ArrayLike[source]

The “+”-branch of the power spectrum of this environment. See the Users Guide on fermionic environments for specifics on the definitions used by QuTiP.

Parameters:
warray_like or float

The frequencies at which to evaluate the power spectrum.

abstract spectral_density(w: float | ArrayLike) float | ArrayLike[source]

The spectral density of this environment. See the Users Guide on fermionic environments for specifics on the definitions used by QuTiP.

Parameters:
warray_like or float

The frequencies at which to evaluate the spectral density.

class LorentzianEnvironment(
T: float,
mu: float,
gamma: float,
W: float,
omega0: float = None,
*,
Nk: int = 10,
tag: Any = None,
)[source]

Bases: FermionicEnvironment

Describes a Lorentzian fermionic environment with the spectral density

\[J(\omega) = \frac{\gamma W^2}{(\omega - \omega_0)^2 + W^2}.\]

(see Eq. 46 in [BoFiN23]).

Parameters:
Tfloat

Environment temperature.

mufloat

Environment chemical potential.

gammafloat

Coupling strength.

Wfloat

The spectral width of the environment.

omega0optional, float (default equal to mu)

The resonance frequency of the environment.

Nkoptional, int, defaults to 10

The number of Pade exponents to be used for the calculation of the correlation functions.

tagoptional, str, tuple or any other object

An identifier (name) for this environment.

approx_by_matsubara(
Nk: int,
tag: Any = None,
) ExponentialFermionicEnvironment[source]

Generates an approximation to this environment by truncating its Matsubara expansion.

Parameters:
Nkint

Number of Matsubara terms to include. In total, the “+” and “-” correlation function branches will include Nk+1 terms each.

tagoptional, str, tuple or any other object

An identifier (name) for the approximated environment. If not provided, a tag will be generated from the tag of this environment.

Returns:
The approximated environment with multi-exponential correlation
function.
approx_by_pade(
Nk: int,
tag: Any = None,
) ExponentialFermionicEnvironment[source]

Generates an approximation to this environment by truncating its Pade expansion.

Parameters:
Nkint

Number of Pade terms to include. In total, the “+” and “-” correlation function branches will include Nk+1 terms each.

tagoptional, str, tuple or any other object

An identifier (name) for the approximated environment. If not provided, a tag will be generated from the tag of this environment.

Returns:
The approximated environment with multi-exponential correlation
function.
correlation_function_minus(
t: float | ArrayLike,
Nk: int = None,
) float | ArrayLike[source]

Calculates the “-“-branch of the two-time auto-correlation function of the Lorentzian environment. The calculation is performed by summing a large number of exponents of the Pade expansion.

Parameters:
tarray_like or float

The time at which to evaluate the correlation function.

Nkint, optional

The number of exponents to use. If not provided, then the value that was provided when the class was instantiated is used.

correlation_function_plus(
t: float | ArrayLike,
Nk: int = None,
) float | ArrayLike[source]

Calculates the “+”-branch of the two-time auto-correlation function of the Lorentzian environment. The calculation is performed by summing a large number of exponents of the Pade expansion.

Parameters:
tarray_like or float

The time at which to evaluate the correlation function.

Nkint, optional

The number of exponents to use. If not provided, then the value that was provided when the class was instantiated is used.

power_spectrum_minus(
w: float | ArrayLike,
) float | ArrayLike[source]

Calculates the “-“-branch of the power spectrum of the Lorentzian environment.

Parameters:
warray_like or float

The frequency at which to evaluate the power spectrum.

power_spectrum_plus(
w: float | ArrayLike,
) float | ArrayLike[source]

Calculates the “+”-branch of the power spectrum of the Lorentzian environment.

Parameters:
warray_like or float

The frequency at which to evaluate the power spectrum.

spectral_density(
w: float | ArrayLike,
) float | ArrayLike[source]

Calculates the Lorentzian spectral density.

Parameters:
warray_like or float

Energy of the mode.

class ExponentialFermionicEnvironment(
ck_plus: ArrayLike = None,
vk_plus: ArrayLike = None,
ck_minus: ArrayLike = None,
vk_minus: ArrayLike = None,
*,
exponents: Sequence[CFExponent] = None,
T: float = None,
mu: float = None,
tag: Any = None,
)[source]

Bases: FermionicEnvironment

Fermionic environment that is specified through an exponential decomposition of its correlation function. The list of coefficients and exponents in the decomposition may either be passed through the four lists ck_plus, vk_plus, ck_minus, vk_minus, or as a list of fermionic CFExponent objects.

Alternative constructors from_plus_exponents and from_minus_exponents are available to compute the “-” exponents automatically from the “+” ones, or vice versa.

Parameters:
ck_pluslist of complex

The coefficients of the expansion terms for the + part of the correlation function. The corresponding frequencies are passed as vk_plus.

vk_pluslist of complex

The frequencies (exponents) of the expansion terms for the + part of the correlation function. The corresponding coefficients are passed as ck_plus.

ck_minuslist of complex

The coefficients of the expansion terms for the - part of the correlation function. The corresponding frequencies are passed as vk_minus.

vk_minuslist of complex

The frequencies (exponents) of the expansion terms for the - part of the correlation function. The corresponding coefficients are passed as ck_minus.

exponentslist of CFExponent

The expansion coefficients and exponents of both parts of the correlation function as CFExponent objects.

T: optional, float

The temperature of the environment.

mu: optional, float

The chemical potential of the environment.

tagoptional, str, tuple or any other object

An identifier (name) for this environment.

correlation_function_minus(
t: float | ArrayLike,
) float | ArrayLike[source]

Computes the “-“-branch of the correlation function represented by this exponential decomposition.

Parameters:
tarray_like or float

The times at which to evaluate the correlation function.

correlation_function_plus(
t: float | ArrayLike,
) float | ArrayLike[source]

Computes the “+”-branch of the correlation function represented by this exponential decomposition.

Parameters:
tarray_like or float

The times at which to evaluate the correlation function.

power_spectrum_minus(
w: float | ArrayLike,
) float | ArrayLike[source]

Calculates the “-“-branch of the power spectrum corresponding to the multi-exponential correlation function.

Parameters:
warray_like or float

The frequency at which to evaluate the power spectrum.

power_spectrum_plus(
w: float | ArrayLike,
) float | ArrayLike[source]

Calculates the “+”-branch of the power spectrum corresponding to the multi-exponential correlation function.

Parameters:
warray_like or float

The frequency at which to evaluate the power spectrum.

spectral_density(
w: float | ArrayLike,
) float | ArrayLike[source]

Computes the spectral density corresponding to the multi-exponential correlation function.

Parameters:
warray_like or float

Energy of the mode.