qml.labs.dla.variational_kak_adj¶
- variational_kak_adj(H, g, dims, adj, verbose=False, opt_kwargs=None, pick_min=False)[source]¶
Variational KaK decomposition of Hermitian
H
using the adjoint representation.Given a Cartan decomposition (
cartan_decomp()
) \(\mathfrak{g} = \mathfrak{k} \oplus \mathfrak{m}\), a Hermitian operator \(H \in \mathfrak{m}\), and a horizontal Cartan subalgebra (horizontal_cartan_subalgebra()
) \(\mathfrak{a} \subset \mathfrak{m}\), this function computes \(a \in \mathfrak{a}\) and \(K_c \in e^{i\mathfrak{k}}\) such that\[H = K_c a K_c^\dagger.\]In particular, \(a = \sum_j c_j a_j\) is decomposed in terms of commuting operators \(a_j \in \mathfrak{a}\). This allows for the immediate decomposition
\[e^{-i t H} = K_c e^{-i t a} K_c^\dagger = K_c \left(\prod_j e^{-i t c_j a_j} \right) K_c^\dagger.\]The result is provided in terms of the adjoint vector representation of \(a \in \mathfrak{a}\) (see
adjvec_to_op()
), i.e. the ordered coefficients \(c_j\) in \(a = \sum_j c_j m_j\) with the basis elements \(m_j \in (\tilde{\mathfrak{m}} \oplus \mathfrak{a})\) and the optimal parameters \(\theta\) such that\[K_c = \prod_{j=|\mathfrak{k}|}^{1} e^{-i \theta_j k_j}\]for the ordered basis of \(\mathfrak{k}\) given by the first \(|\mathfrak{k}|\) elements of
g
. Note that we define \(K_c\) mathematically with the descending order of basis elements \(k_j \in \mathfrak{k}\) such that the resulting circuit has the canonical ascending order. In particular, a PennyLane quantum function that describes the circuit given the optimal parameterstheta_opt
and the basisk
containing the operators, is given by the following.def Kc(theta_opt: Iterable[float], k: Iterable[Operator]): assert len(theta_opt) == len(k) for theta_j, k_j in zip(theta_opt, k): qml.exp(-1j * theta_j * k_j)
Internally, this function performs a modified version of 2104.00728, in particular minimizing the cost function
\[f(\theta) = \langle H, K(\theta) e^{-i \sum_{j=1}^{|\mathfrak{a}|} \pi^j a_j} K(\theta)^\dagger \rangle,\]see eq. (6) therein and our demo for more details. Instead of relying on having Pauli words, we use the adjoint representation for a more general evaluation of the cost function. The rest is the same.
- Parameters:
H (Union[Operator, PauliSentence, np.ndarray]) – Hamiltonian to decompose
g (List[Union[Operator, PauliSentence, np.ndarray]]) – DLA of the Hamiltonian
dims (Tuple[int]) – Tuple of dimensions
(dim_k, dim_mtilde, dim_a)
of Cartan decomposition \(\mathfrak{g} = \mathfrak{k} \oplus (\tilde{\mathfrak{m}} \oplus \mathfrak{a})\)adj (np.ndarray) – Adjoint representation of dimension
(dim_g, dim_g, dim_g)
, with the implicit ordering(k, mtilde, a)
.verbose (bool) – Plot the optimization. Requires matplotlib to be installed (
pip install matplotlib
)opt_kwargs (dict) – Keyword arguments for the optimization like initial starting values for \(\theta\) of dimension
(dim_k,)
, given astheta0
. Also includesn_epochs
,lr
,b1
,b2
,verbose
,interrupt_tol
, seerun_opt()
pick_min (bool) – Whether to pick the parameter set with lowest cost function value during the optimization as optimal parameters. Otherwise picks the last parameter set.
- Returns:
(adjvec_a, theta_opt)
: The adjoint vector representationadjvec_a
of dimension(dim_mtilde + dim_a,)
, with respect to the basis of \(\mathfrak{m} = \tilde{\mathfrak{m}} + \mathfrak{a}\) of the CSA element \(a \in \mathfrak{a}\) s.t. \(H = K a K^\dagger\). For a successful optimization, the entries corresponding to \(\tilde{\mathfrak{m}}\) should be close to zero. The second return value,theta_opt
, are the optimal coefficients \(\theta\) of the decomposition \(K = \prod_{j=|\mathfrak{k}|}^{1} e^{-i \theta_j k_j}\) for the basis \(k_j \in \mathfrak{k}\).- Return type:
Tuple(np.ndarray, np.ndarray)
Example
Let us perform a KaK decomposition for the transverse field Ising model Hamiltonian, exemplarily for \(n=3\) qubits on a chain. We start with some boilerplate code to perform a Cartan decomposition using the
concurrence_involution()
, which places the Hamiltonian in the horizontal subspace \(\mathfrak{m}\). From this we re-order \(\mathfrak{g} = \mathfrak{k} + \mathfrak{m}\) and finally compute ahorizontal_cartan_subalgebra()
\(\mathfrak{a}\) in \(\mathfrak{m} = \tilde{\mathfrak{m}} \oplus \mathfrak{a}\).import pennylane as qml import numpy as np import jax.numpy as jnp import jax from pennylane import X, Z from pennylane.liealg import ( cartan_decomp, horizontal_cartan_subalgebra, check_cartan_decomp, concurrence_involution, adjvec_to_op, ) from pennylane.labs.dla import ( validate_kak, variational_kak_adj, ) n = 3 gens = [X(i) @ X(i + 1) for i in range(n - 1)] gens += [Z(i) for i in range(n)] H = qml.sum(*gens) g = qml.lie_closure(gens) g = [op.pauli_rep for op in g] involution = concurrence_involution assert not involution(H) k, m = cartan_decomp(g, involution=involution) assert check_cartan_decomp(k, m) g = k + m adj = qml.structure_constants(g) g, k, mtilde, a, adj = horizontal_cartan_subalgebra(g, k, m, adj, tol=1e-14, start_idx=0)
Due to the canonical ordering of all constituents, it suffices to tell
variational_kak_adj
the dimensions ofdims = (len(k), len(mtilde), len(a))
, alongside the HamiltonianH
, the Lie algebrag
and its adjoint representationadj
. Internally, the function is performing a variational optimization to find a local extremum of a suitably constructed loss function that finds as its extremum the decomposition\[K_c = \prod_{j=1}^{|\mathfrak{k}|} e^{-i \theta_j k_j}\]in form of the optimal parameters \(\{\theta_j\}\) for the respective \(k_j \in \mathfrak{k}\). The resulting \(K\) then informs the CSA element
a
of the KaK decomposition via \(a = K_c H K_c^\dagger\). This is detailed in 2104.00728.>>> dims = (len(k), len(mtilde), len(a)) >>> adjvec_a, theta_opt = variational_kak_adj(H, g, dims, adj, opt_kwargs={"n_epochs": 3000})
As a result, we are provided with the adjoint vector representation of the CSA element \(a \in \mathfrak{a}\) with respect to the basis
mtilde+a
and the optimal parameters of dimension \(|\mathfrak{k}|\)Let us perform some sanity checks to better understand the resulting outputs. We can turn that element back to an operator using
adjvec_to_op()
and from that to a matrix for which we can check Hermiticity.m = mtilde + a [a_op] = adjvec_to_op([adjvec_a], m) a_m = qml.matrix(a_op, wire_order=range(n)) assert np.allclose(a_m, a_m.conj().T)
Let us now confirm that we get back the original Hamiltonian from the resulting \(K_c\) and \(a\). In particular, we want to confirm \(H = K_c a K_c^\dagger\) for \(K_c = \prod_{j=1}^{|\mathfrak{k}|} e^{-i \theta_j k_j}\).
assert len(theta_opt) == len(k) def Kc(theta_opt): for th, op in zip(theta_opt, k): qml.exp(-1j * th * op.operation()) Kc_m = qml.matrix(Kc, wire_order=range(n))(theta_opt) # check Unitary property of Kc assert np.allclose(Kc_m.conj().T @ Kc_m, np.eye(2**n)) H_reconstructed = Kc_m @ a_m @ Kc_m.conj().T H_m = qml.matrix(H, wire_order=range(len(H.wires))) # check Hermitian property of reconstructed Hamiltonian assert np.allclose( H_reconstructed, H_reconstructed.conj().T ) # confirm reconstruction was successful to some given numerical tolerance assert np.allclose(H_m, H_reconstructed, atol=1e-6)
Instead of performing these checks by hand, we can use the helper function
validate_kak()
.>>> assert validate_kak(H, g, k, (adjvec_a, theta_opt), n, 1e-6)