rydiqule.solvers.solve_steady_state¶
- rydiqule.solvers.solve_steady_state(sensor: Sensor, doppler: bool = False, doppler_mesh_method: UniformMethod | IsoPopMethod | SplitMethod | DirectMethod | None = None, sum_doppler: bool = True, weight_doppler: bool = True, n_slices: int | None = None) Solution [source]¶
Finds the steady state solution for a system characterized by a sensor.
If insuffucent system memory is available to solve the system in a single call, system is broken into “slices” of manageable memory footprint which are solved indivudually. This slicing behavior does not affect the result. Can be performed with or without doppler averging.
- Parameters:
sensor (
Sensor
) – The sensor for which the solution will be calculated.doppler (bool, optional) – Whether to calculate the solution for a doppler-broadened gas. If
True
, only uses dopper brodening defined bykvec
parameters for couplings in thesensoe
, so setting thisTrue
withoutkvec
definitions will have no effect. Default isFalse
.(dict (doppler_mesh_method) – If not
None
, should be a dictionary of meshing parameters to be passed todoppler_classes()
. Seedoppler_classes()
for more information on supported methods and arguments. IfNone, uses the default doppler meshing. Default is `None
.optional) – If not
None
, should be a dictionary of meshing parameters to be passed todoppler_classes()
. Seedoppler_classes()
for more information on supported methods and arguments. IfNone, uses the default doppler meshing. Default is `None
.sum_doppler (bool) – Whether to average over doppler classes after the solve is complete. Setting to
False
will not perform the sum, allowing viewing of the weighted results of the solve for each doppler class. In this case, an axis will be prepended to the solution for each axis along which doppler broadening is computed. Ignored ifdoppler=False
. Default isTrue
.weight_doppler (bool) – Whether to apply weights to doppler solution to perform averaging. If
False
, will not apply weights or perform a doppler_average, regardless of the value ofsum_doppler
. Changing from default intended only for internal use. Ignored ifdoppler=False
orsum_doppler=False
. Default isTrue
.n_slices (int or None, optional) – How many sets of equations to break the full equations into. The actual number of slices will be the largest between this value and the minumum number of slices to solve the system without a memory error. If
None
, uses the minimum number of slices to solve the system without a memory error. Detailed information about slicing behavior can be found inmatrix_slice()
. Default isNone
.
Notes
Note
If decoherence values are not sufficiently populated in the sensor, the resulting equations may be singular, resulting in an error in
numpy.linalg
. This error is not caught for flexibility, but is likely the culprit fornumpy.linalg
errors encountered in steady-state solves.Note
The solution produced by this function will be expressed using rydiqule’s convention of converting a density matrix into the real basis and removing the ground state to improve numerical stability.
Note
If the sensor contains couplings with
time_dependence
, this solver will add those couplings at their \(t=0\) value to the steady-state hamiltonian to solve.- Returns:
An object contining the solution and related information.
- Return type:
Examples
A basic solve for a 3-level system would have a “density matrix” solution of size 8 (3^2-1)
>>> s = rq.Sensor(3) >>> s.add_coupling((0,1), detuning = 1, rabi_frequency=1) >>> s.add_coupling((1,2), detuning = 2, rabi_frequency=2) >>> s.add_transit_broadening(0.1) >>> sol = rq.solve_steady_state(s) >>> print(type(sol)) <class 'rydiqule.sensor_solution.Solution'> >>> print(type(sol.rho)) <class 'numpy.ndarray'> >>> print(sol.rho.shape) (8,)
Defining an array-like parameter will automatically calculate the density matrix solution for every value. Here we use 11 values, resulting in 11 density matrices. The
axis_labels
attribute of the solution can clarify which axes are which.>>> s = rq.Sensor(3) >>> det = np.linspace(-1,1,11) >>> s.add_coupling((0,1), detuning = det, rabi_frequency=1) >>> s.add_coupling((1,2), detuning = 2, rabi_frequency=2) >>> s.add_transit_broadening(0.1) >>> sol = rq.solve_steady_state(s) >>> print(sol.rho.shape) (11, 8) >>> print(sol.axis_labels) ['(0,1)_detuning', 'density_matrix']
>>> s = rq.Sensor(3) >>> det = np.linspace(-1,1,11) >>> s.add_coupling((0,1), detuning = det, rabi_frequency=1) >>> s.add_coupling((1,2), detuning = det, rabi_frequency=2) >>> s.add_transit_broadening(0.1) >>> sol = rq.solve_steady_state(s) >>> print(sol.rho.shape) (11, 11, 8) >>> print(sol.axis_labels) ['(0,1)_detuning', '(1,2)_detuning', 'density_matrix']
If the solve uses doppler broadening, but not averaging for doppler is specified, there will be a solution axis corresponding to doppler classes.
>>> s = rq.Sensor(3, vP=1) >>> det = np.linspace(-1,1,11) >>> s.add_coupling((0,1), detuning = det, rabi_frequency=1) >>> s.add_coupling((1,2), detuning = 2, rabi_frequency=2, kvec=(4,0,0)) >>> s.add_transit_broadening(0.1) >>> sol = rq.solve_steady_state(s, doppler=True, sum_doppler=False) >>> print(sol.rho.shape) (561, 11, 8) >>> print(sol.axis_labels) ['doppler_0', '(0,1)_detuning', 'density_matrix']