qml.liealg.horizontal_cartan_subalgebra¶
- horizontal_cartan_subalgebra(k, m, adj=None, start_idx=0, tol=1e-10, verbose=0, return_adjvec=False, is_orthogonal=True)[source]¶
Compute a Cartan subalgebra (CSA) a⊆m.
A non-unique CSA is a maximal Abelian subalgebra in the horizontal subspace m of a Cartan decomposition. Note that this is sometimes called a horizontal CSA, and is different from the other definitions of a CSA.
The final decomposition yields
g=k⊕(˜m⊕a),where a) is the CSA and ˜m is the remainder of the horizontal subspace m.
See also
cartan_decomp()
,structure_constants()
, The KAK decomposition in theory (demo), The KAK decomposition in practice (demo).- Parameters
k (List[Union[PauliSentence, TensorLike]]) – Vertical space k from Cartan decomposition g=k⊕m.
m (List[Union[PauliSentence, TensorLike]]) – Horizontal space m from Cartan decomposition g=k⊕m.
adj (Array) – The |g|×|g|×|g| dimensional adjoint representation of g. When
None
is provided,structure_constants()
is used internally by default to compute the adjoint representation.start_idx (bool) – Indicates from which element in
m
the CSA computation starts.tol (float) – Numerical tolerance for linear independence check.
verbose (bool) – Whether or not to output progress during computation.
return_adjvec (bool) – Determine the output format. If
False
, returns operators in their original input format (matrices orPauliSentence
). IfTrue
, returns the spaces as adjoint representation vectors (see :func`~op_to_adjvec` and ~adjvec_to_op).is_orthogonal (bool) – Whether the basis elements are all orthogonal, both within and between
g
,k
andm
.
- Returns
A tuple of adjoint vector representations
(newg, k, mtilde, a, new_adj)
, corresponding to g, k, ˜m, a and the new adjoint representation. The dimensions are(|g|, |g|)
,(|k|, |g|)
,(|mtilde|, |g|)
,(|a|, |g|)
and(|g|, |g|, |g|)
, respectively.- Return type
Tuple(TensorLike, TensorLike, TensorLike, TensorLike, TensorLike)
Example
A quick example computing a Cartan subalgebra of su(4) using the Cartan involution
even_odd_involution()
.>>> import pennylane as qml >>> g = list(qml.pauli.pauli_group(2)) # u(4) >>> g = g[1:] # remove identity -> su(4) >>> g = [op.pauli_rep for op in g] # optional; turn to PauliSentence for convenience >>> k, m = qml.liealg.cartan_decomp(g, qml.liealg.even_odd_involution) >>> g = k + m # re-order g to separate k and m >>> newg, k, mtilde, a, new_adj = qml.liealg.horizontal_cartan_subalgebra(k, m) >>> newg == k + mtilde + a True >>> a [-1.0 * Z(0) @ Z(1), 1.0 * Y(0) @ Y(1), -1.0 * X(0) @ X(1)]
We can confirm that these all commute with each other, as the CSA is Abelian (i.e., all operators commute).
>>> qml.liealg.check_abelian(a) True
We can opt-in to return what we call adjoint vectors of dimension |g|, where each component corresponds to an entry in (the ordered)
g
. The adjoint vectors for the Cartan subalgebra are innp_a
.from pennylane.liealg import horizontal_cartan_subalgebra np_newg, np_k, np_mtilde, np_a, new_adj = horizontal_cartan_subalgebra(k, m, return_adjvec=True)
We can reconstruct an operator by computing ˆOv=∑ivigi for an adjoint vector v and gi∈g.
>>> v = np_a[0] >>> op = sum(v_i * g_i for v_i, g_i in zip(v, g)) >>> op.simplify() >>> op -1.0 * Z(0) @ Z(1)
For convenience, we provide a helper function
adjvec_to_op()
for conversion of the returned collections of adjoint vectors.>>> a = qml.liealg.adjvec_to_op(np_a, g) >>> a [-1.0 * Z(0) @ Z(1), 1.0 * Y(0) @ Y(1), -1.0 * X(0) @ X(1)]
Usage Details
Let us walk through an example of computing the Cartan subalgebra. The basis for computing the Cartan subalgebra is having the Lie algebra g, a Cartan decomposition g=k⊕m and its adjoint representation.
We start by computing these ingredients using
cartan_decomp()
andstructure_constants()
. As an example, we take the Lie algebra of the Heisenberg model with generators {XiXi+1,YiYi+1,ZiZi+1}.>>> from pennylane.liealg import cartan_decomp >>> from pennylane import X, Y, Z >>> n = 3 >>> gens = [X(i) @ X(i+1) for i in range(n-1)] >>> gens += [Y(i) @ Y(i+1) for i in range(n-1)] >>> gens += [Z(i) @ Z(i+1) for i in range(n-1)] >>> g = qml.lie_closure(gens, matrix=True)
Taking the Heisenberg Lie algebra, we can perform the Cartan decomposition. We take the
even_odd_involution()
as a valid Cartan involution. The resulting vertical and horizontal subspaces k and m need to fulfill the commutation relations [k,k]⊆k, [k,m]⊆m and [m,m]⊆k, which we can check using the helper functioncheck_cartan_decomp()
.>>> from pennylane.liealg import even_odd_involution, check_cartan_decomp >>> k, m = cartan_decomp(g, even_odd_involution) >>> check_cartan_decomp(k, m) # check commutation relations of Cartan decomposition True
Our life is easier when we use a canonical ordering of the operators. This is why we re-define
g
with the new ordering in terms of operators ink
first, and then all remaining operators fromm
.>>> import numpy as np >>> g = np.vstack([k, m]) # re-order g to separate k and m operators >>> adj = qml.structure_constants(g, matrix=True) # compute adjoint representation of g
Finally, we can compute a Cartan subalgebra a, a maximal Abelian subalgebra of m.
>>> newg, k, mtilde, a, new_adj = horizontal_cartan_subalgebra(k, m, adj, start_idx=3)
The new DLA
newg
is just the concatenation ofk
,mtilde
,a
. Each component is returned in the original input format. Here we obtain collections of 8×8 matrices (numpy
arrays), as this is what we started from.>>> newg.shape, k.shape, mtilde.shape, a.shape, new_adj.shape ((15, 8, 8), (6, 8, 8), (6, 8, 8), (3, 8, 8), (15, 15, 15))
We can also let the function return what we call adjoint representation vectors.
>>> kwargs = {"start_idx": 3, "return_adjvec": True} >>> np_newg, np_k, np_mtilde, np_a, new_adj = horizontal_cartan_subalgebra(k, m, adj, **kwargs) >>> np_newg.shape, np_k.shape, np_mtilde.shape, np_a.shape, new_adj.shape ((15, 15), (6, 15), (6, 15), (3, 15), (15, 15, 15))
These are dense vector representations of dimension |g|, in which each entry corresponds to the respective operator in g. Given an adjoint representation vector v, we can reconstruct the respective operator by simply computing ˆOv=∑ivigi, where gi∈g (hence the need for a canonical ordering).
We provide a convenience function
adjvec_to_op()
that works with bothg
represented as dense matrices or PL operators. Because we used dense matrices in this example, we transform the operators back to PennyLane operators usingpauli_decompose()
.>>> from pennylane.liealg import adjvec_to_op >>> a = adjvec_to_op(np_a, g) >>> h_op = [qml.pauli_decompose(op).pauli_rep for op in a] >>> h_op [-1.0 * Y(1) @ Y(2), -1.0 * Z(1) @ Z(2), 1.0 * X(1) @ X(2)]
In that case we chose a Cartan subalgebra from which we can readily see that it is commuting, but we also provide a small helper function to check that.
>>> from pennylane.liealg import check_abelian >>> check_abelian(h_op) True
Last but not least, the adjoint representation
new_adj
is updated to represent the new basis and its ordering ofg
.