qdyn.model module¶
Abstraction for complete physical models that can be written to a runfolder as a config file and dependent data
Summary¶
Classes:
Model for a system well-described in the energy basis. |
Functions:
Return a random N x N density matrix. |
Reference¶
-
qdyn.model.
random_density_matrix
(N)[source]¶ Return a random N x N density matrix.
>>> rho = random_density_matrix(10)
The resulting density matrix is normalized, positive semidefinite, and Hermitian
>>> assert( abs(np.trace(rho) - 1.0) <= 1.0e-14 ) >>> assert( np.min(np.linalg.eigvals(rho).real) > 0.0 ) >>> assert( np.max(np.abs(rho.conjugate().transpose() - rho)) <= 1.0e-14 )
-
class
qdyn.model.
LevelModel
[source]¶ Bases:
object
Model for a system well-described in the energy basis. That is, all operators are (sparse) matrices, and all states are simple vectors
-
t0
¶ Initial time.
- Type
-
T
¶ Final time.
- Type
-
use_mcwf
¶ Propagate using the Monte-Carlo Wave Function (quantum jump) method
-
construct_mcwf_ham
¶ When using the MCWF method, the propagation must use an “effective” Hamiltonian that includes a non-Hermitian decay term. This term is constructed from the Lindblad operators. If
use_mcwf=True
andconstruct_mcwf_ham=False
, it is the user’s responsibility to ensure that ham is the proper effective Hamiltonian. The construct_mcwf flag determines the presence of add_to_H_jump config file parameter for each Lindblad operator- Type
-
user_data
¶ Key-value pairs that should that describe user-defined data. These will go in the
user_strings
,user_reals
,user_logicals
, oruser_ints
section of the config file, depending on the type of the value
After instantiation, the attributes t0, T, nt, prop_method, use-mcwf, and construct_mcwf_ham are all set via
set_propagation()
. Operators and pulses are added to the system throughadd_ham()
,add_observable()
, andadd_lindblad_op()
. States are added throughadd_state()
. Both the general OCT settings (OCT section in the QDYN config file) and OCT-related settings for each control pulse are controlled throughset_oct()
. After the model has been constructed, a config file and all dependent data input files for the operators, pulses, and states can be written viawrite_to_runfolder()
.-
observables
(label=None, with_attribs=False)[source]¶ Return list of all observables with the matching label (or all labels if label is ‘*’). If with_attribs is True, the result is a list of tuples
(observable, attributes
) whereattributes
is a dictionary of config file attributes
-
lindblad_ops
(label=None, with_attribs=False)[source]¶ Return list of all Lindblad operators with the matching label (or all labels if label is ‘*’). If with_attribs is True, the result is a list of tuples
(operator, attributes
) whereattributes
is a dictionary of config file attributes
-
ham
(label=None, with_attribs=False)[source]¶ Return list of all Hamiltonian operators with the matching label (or all labels if label is ‘*’). If with_attribs is True, the result is a list of tuples
(ham, attributes
) whereattributes
is a dictionary of config file attributes
-
pulses
(label=None)[source]¶ Return a list of all known pulses with the matching label (or all labels if label is ‘*’).
-
add_ham
(H, pulse=None, op_unit=None, sparsity_model=None, op_type=None, label=None, **kwargs)[source]¶ Add a term to the system Hamiltonian. If called repeatedly, the total Hamiltonian is the sum of all added terms.
- Parameters
H – Hamiltonian matrix. Can be a numpy matrix or array, scipy sparse matrix, or qutip.Qobj
pulse – if not None, H will couple to pulse. Can be an instance of
AnalyticalPulse
orPulse
.sparsity_model (None or str) – sparsity model that QDYN should use to encode the data in H. If None, will be determined automatically
op_type (None or str) – the value for
op_type
in the config file that should be used for H. This determines how exactly H couples to pulse. Common values are ‘dipole’ (linear coupling) and ‘dstark’ (quadratic “Stark shift” coupling)label (str or None) – Multiple Hamiltonians may be defined in the same config file if they are differentiated by label. The default label is the empty string
kwargs – All other keyword arguments set options for H in the config file (e.g. specrad_method, filename)
Note
It is recommended to use
AnalyticalPulse
to express time-dependency, as this is independent of a specific time grid. Instances ofPulse
must have a time grid that exactly matches the time grid specified viaset_propagation()
.
-
add_observable
(O, outfile, exp_unit, time_unit, col_label, square=None, exp_surf=None, is_real=None, in_lab_frame=False, op_unit=None, label=None, pulse=None, from_time_index=None, to_time_index=None, step=None)[source]¶ Add an observable
- Parameters
O (numpy.matrix, str) – Observable to add. Must be a matrix, or one of “ham”, “norm”, “pop”
outfile (str) – Name of output file to which to write expectation values of O
exp_unit (str) – The unit in which to write the expectation value in outfile, as well as the default value for op_unit
time_unit (str) – The unit in which to write the time grid in outfile
col_label (str) – The label for the column in outfile containing the expectation value of O
square (str or None) – If not None, label for the column in outfile containing the expectation value for the square of O
exp_surf (int or None) – The surface number the expectation value; only if O is a string
is_real (bool or None) – Whether or not the expectation value is real. If not given, this is determined automatically
in_lab_frame (bool) – If True, indicates that the observable is defined in the lab frame. When expectation values are calculated, this should not be done with states in the rotating frame.
op_unit (str or None) – The unit in which the entries of O are given. By default, this is the same as exp_unit.
label (str or None) – Observables associated with different Hamiltonians may be defined in the same config file if they are differentiated by label. The default label is the empty string.
pulse – If not None, a pulse for the observable to couple to (see
add_ham()
)from_time_index (int or None) – Index on total time grid from where to start calculating/writing expectation values. Negative values count from the end of the time grid. If None, equivalent to 1.
to_time_index (int or None) – Index on total time grid where to stop calculating/writing expectation values. Negative values count from the end of the time grid. If None, equivalent to -1
step (int or None) – Step width for where to calculate/write expectation values between from_time_index and to_time_index. If None, equivalent to 1
-
add_lindblad_op
(L, op_unit=None, sparsity_model=None, add_to_H_jump=None, label=None, pulse=None, **kwargs)[source]¶ Add Lindblad operator.
- Parameters
L (tuple, numpy.matrix) – Lindblad operator to add. Must be a matrix or a tuple
(matrix, pulse)
, cf. the ham attribute.op_unit (None or str) – Unit of the values in L, e.g.
sqrt_GHz
(Lindblad operators are in units square-root-of-energy)sparsity_model (None or str) – sparsity model that QDYN should use to encode the data in L. If None, will be determined automatically
add_to_H_jump (None, str) – sparsity model to be set for add_to_H_jump, i.e. for the decay term in the effective MCWF Hamiltonian. If None, will be de determined automatically.
label (str or None) – Lindblad operators associated with different Hamiltonians may be defined in the same config file if they are differentiated by label. The default label is the empty string.
pulse – If not None, a pulse for the Lindblad operator to couple to (see
add_ham()
)
All other keyword arguments set options for L in the config file.
-
set_dissipator
(D, op_unit=None, sparsity_model=None, label=None, pulse=None, **kwargs)[source]¶ Set a dissipation superoperator in the config file
- Parameters
D (tuple, numpy.matrix) – Dissiption superoperoperator to add. Must be a matrix or a tuple
(matrix, pulse)
op_unit (None or str) – Unit of the value in D. e.g.
GHz
(dissipators are in units of energy)sparsity_model (None or str) – sparsity model that QDYN should use to encode the data in D. If None, will be determined automatically
label (str or None) – Dissipators associated with different Hamiltonians may be defined in the same config file if they are differentiated by label. The default label is the empty string.
pulse – If not None, a pulse for the Dissipator to couple to
All other keyword arguments set options for the dissipator in the config file.
-
set_propagation
(T, nt, time_unit, t0=0.0, prop_method='newton', use_mcwf=False, mcwf_order=2, construct_mcwf_ham=True, label=None, initial_state=None)[source]¶ Set the time grid and other propagation-specific settings
- Parameters
T (float) – final time
nt (int) – number of points in the time grid
time_unit (str) – physical unit of T, t0
t0 (float) – initial time
prop_method (str) – method to be used for propagation
use_mcwf (bool) – If True, prepare for Monte-Carlo wave function propagation
mcwf_order (int) – Order for MCWF, must be 1 or 2
construct_mcwf_ham (bool) – When using MCWF (use_mcwf=True), by default an additional inhomogeneous decay term is added to the Hamiltonian, based on the Lindblad operators. By passing construct_mcwf_ham=False, this does not happen. It is the user’s responsibility then to ensure the Hamiltonian in the model is the correct “effective” Hamiltonian for a MCWF propagation.
initial_state (numpy.ndarray or None) – Initial wave function
Notes
When setting up an MCWF propagation, using the mcwf_order=2 is usually the right thing to do. In some cases of strong dissipation, it may be numerically more efficient to use a first-order MCWF method, where in each time interval dt between two time grids there is at most one quantum jump, and the jumps takes place over the entire duration dt. The time steps must be very small accordingly. In contrast, for mcwf_order=2, all jumps are instantaneous, and there can be multiple jumps per time step, but the numerical overhead is larger.
-
set_oct
(method, J_T_conv, max_ram_mb, **kwargs)[source]¶ Set config file data and pulse properties for optimal control
- Parameters
method (str) – Optimization method. Allowed values are ‘krotovpk’, ‘krotov2’, ‘krotovbwr’, and ‘lbfgs’
J_T_conv (float) – The value of the final time functional
max_ram_mb (int) – The amount of memory that is available for storing propagated states. If this is not suffient to hold all the states required to calculate the gradient, a “segmented” storage scheme will be used that caches the states to disk, using up to max_disk_mb hard drive storage.
All other keyword arguments directly specify keys and values for the OCT config section. Allowed keys are iter_start, iter_stop, max_disk_mb, lbfgs_memory, linesearch, grad_order, iter_dat, tau_dat, params_file, krotov2_conv_dat, ABC_dat, delta_J_conv, delta_J_T_conv, A, B, C, dynamic_sigma, dynamic_lambda_a, strict_convergence, limit_pulses, sigma_form, max_seconds, lambda_b, keep_pulses, re_init_prop, continue storage_folder, bwr_nint, bwr_base, g_a, see the QDYN Fortran library documentation for details.
Note: continue is a python keyword, which is why the key for continue is can also be set with continue_
- Raises
KeyError – If the settings for a pulse contain invalid or missing keys
-