weylchamber.local_invariants module¶
Summary¶
Functions:
J_T_LI |
Calculate value of the local-invariants functional |
closest_LI |
Find the closest gate that has the given Weyl chamber coordinates |
g1g2g3 |
Calculate local invariants \((g_1, g_3, g_3)\) |
g1g2g3_from_c1c2c3 |
Calculate local invariants from the Weyl chamber coordinates |
make_LI_krotov_chi_constructor |
Return a constructor for the χ’s in an LI optimization. |
__all__
: J_T_LI
, closest_LI
, g1g2g3
, g1g2g3_from_c1c2c3
, make_LI_krotov_chi_constructor
Reference¶
-
weylchamber.local_invariants.
g1g2g3
(U, ndigits=8)[source]¶ Calculate local invariants \((g_1, g_3, g_3)\)
Given a two-qubit gate, calculate local invariants \((g_1, g_2, g_3)\). U must be in the canonical basis. For numerical stability, the resulting values are rounded to the given precision, cf. the ndigits parameter of the built-in
round()
function.>>> print("%.2f %.2f %.2f" % g1g2g3(qutip.gates.cnot())) 0.00 0.00 1.00
Return type: Tuple
[float
,float
,float
]
-
weylchamber.local_invariants.
g1g2g3_from_c1c2c3
(c1, c2, c3, ndigits=8)[source]¶ Calculate local invariants from the Weyl chamber coordinates
Calculate the local invariants \((g_1, g_2, g_3)\) from the Weyl chamber coordinates \((c_1, c_2, c_3)\), in units of π. The result is rounded to the given precision, in order to enhance numerical stability (cf. ndigits parameter of the built-in
round()
function)Example
>>> CNOT = qutip.gates.cnot() >>> print("%.2f %.2f %.2f" % g1g2g3_from_c1c2c3(*c1c2c3(CNOT))) 0.00 0.00 1.00
Return type: Tuple
[float
,float
,float
]
-
weylchamber.local_invariants.
J_T_LI
(O, U, form='g')[source]¶ Calculate value of the local-invariants functional
Parameters:
-
weylchamber.local_invariants.
closest_LI
(U, c1, c2, c3, method='leastsq', limit=1e-06)[source]¶ Find the closest gate that has the given Weyl chamber coordinates
The c1, c2, c3 are given in units of π
-
weylchamber.local_invariants.
make_LI_krotov_chi_constructor
(gate, canonical_basis, unitarity_weight=0)[source]¶ Return a constructor for the χ’s in an LI optimization.
Return a chi_constructor that determines the boundary condition of the backwards propagation in an optimization towards the local equivalence class of gate in Krotov’s method, based on the foward-propagtion of the Bell states.
Parameters: - gate (qutip.Qobj) – A 4×4 quantum gate, in the canonical_basis.
- canonical_basis (list[qutip.Qobj]) – A list of four basis states that define the canonical basis \(\ket{00}\), \(\ket{01}\), \(\ket{10}\), and \(\ket{11}\) of the logical subspace.
- unitarity_weight (float) – A weight in [0, 1] that determines how much emphasis is placed on maintaining population in the logical subspace.
Returns: a function
chi_constructor(fw_states_T, *args)
that receive the result of a foward propagation of the Bell states (obtained from canonical_basis viaweylchamber.gates.bell_basis()
), and returns a list of statex \(\ket{\chi}\) that are the boundary condition for the backward propagation in Krotov’s method. Positional arguments beyond fw_states_T are ignored.Return type: callable