diff --git a/python/ffsim/__init__.py b/python/ffsim/__init__.py index 2686a6c5e..fd9690c53 100644 --- a/python/ffsim/__init__.py +++ b/python/ffsim/__init__.py @@ -15,6 +15,7 @@ from ffsim.gates import ( apply_diag_coulomb_evolution, apply_fsim_gate, + apply_fswap_gate, apply_givens_rotation, apply_hop_gate, apply_num_interaction, @@ -132,6 +133,7 @@ "addresses_to_strings", "apply_diag_coulomb_evolution", "apply_fsim_gate", + "apply_fswap_gate", "apply_givens_rotation", "apply_hop_gate", "apply_num_interaction", diff --git a/python/ffsim/gates/__init__.py b/python/ffsim/gates/__init__.py index 6d530796b..993f0eecf 100644 --- a/python/ffsim/gates/__init__.py +++ b/python/ffsim/gates/__init__.py @@ -12,6 +12,7 @@ from ffsim.gates.basic_gates import ( apply_fsim_gate, + apply_fswap_gate, apply_givens_rotation, apply_hop_gate, apply_num_interaction, @@ -27,6 +28,7 @@ __all__ = [ "apply_diag_coulomb_evolution", "apply_fsim_gate", + "apply_fswap_gate", "apply_givens_rotation", "apply_hop_gate", "apply_num_interaction", diff --git a/python/ffsim/gates/basic_gates.py b/python/ffsim/gates/basic_gates.py index ea9fc6b16..a9b101baa 100644 --- a/python/ffsim/gates/basic_gates.py +++ b/python/ffsim/gates/basic_gates.py @@ -570,3 +570,60 @@ def apply_fsim_gate( vec, -phi, target_orbs, norb=norb, nelec=nelec, spin=spin, copy=False ) return vec + + +def apply_fswap_gate( + vec: np.ndarray, + target_orbs: tuple[int, int], + norb: int, + nelec: int | tuple[int, int], + spin: Spin = Spin.ALPHA_AND_BETA, + *, + copy: bool = True, +) -> np.ndarray: + r"""Apply an FSWAP gate. + + The FSWAP gate swaps two orbitals. It is represented by the operator + + .. math:: + \text{FSWAP}(p, q) = + 1 + a^\dagger_p a_q + a^\dagger_q a_p - a_p^\dagger a_p - a_q^\dagger a_q + + Under the Jordan-Wigner transform, this gate has the following matrix when applied + to neighboring qubits: + + .. math:: + + \begin{pmatrix} + 1 & 0 & 0 & 0 \\ + 0 & 0 & 1 & 0 \\ + 0 & 1 & 0 & 0 \\ + 0 & 0 & 0 & -1 \\ + \end{pmatrix} + + Args: + vec: The state vector to be transformed. + target_orbs: The orbitals (p, q) to swap. + norb: The number of spatial orbitals. + nelec: Either a single integer representing the number of fermions for a + spinless system, or a pair of integers storing the numbers of spin alpha + and spin beta fermions. + spin: Choice of spin sector(s) to act on. + copy: Whether to copy the vector before operating on it. + + - If `copy=True` then this function always returns a newly allocated + vector and the original vector is left untouched. + - If `copy=False` then this function may still return a newly allocated + vector, but the original vector may have its data overwritten. + It is also possible that the original vector is returned, + modified in-place. + """ + if copy: + vec = vec.copy() + mat = np.eye(norb, dtype=complex) + mat[np.ix_(target_orbs, target_orbs)] = [[0, 1], [1, 0]] + if isinstance(nelec, int): + return apply_orbital_rotation(vec, mat, norb=norb, nelec=nelec, copy=False) + return apply_orbital_rotation( + vec, pair_for_spin(mat, spin=spin), norb=norb, nelec=nelec, copy=False + ) diff --git a/tests/python/gates/basic_gates_test.py b/tests/python/gates/basic_gates_test.py index f5b55bbd9..5260f7472 100644 --- a/tests/python/gates/basic_gates_test.py +++ b/tests/python/gates/basic_gates_test.py @@ -595,3 +595,57 @@ def mat(theta: float) -> np.ndarray: phase_11=phase_11, norb=norb, ) + + +@pytest.mark.parametrize("norb, spin", ffsim.testing.generate_norb_spin(range(4))) +def test_apply_fswap_gate_matrix_spinful(norb: int, spin: ffsim.Spin): + """Test applying fSWAP gate matrix.""" + + mat01 = np.array([[0, 1], [1, 0]]) + + phase_00 = 1 + phase_11 = -1 + + for i, j in itertools.combinations(range(norb), 2): + for target_orbs in [(i, j), (j, i)]: + assert_has_two_orbital_matrix( + lambda vec, norb, nelec: ffsim.apply_fswap_gate( + vec, + target_orbs=target_orbs, + norb=norb, + nelec=nelec, + spin=spin, + ), + target_orbs=target_orbs, + mat=mat01, + phase_00=phase_00, + phase_11=phase_11, + norb=norb, + spin=spin, + ) + + +@pytest.mark.parametrize("norb", range(4)) +def test_apply_fswap_gate_matrix_spinless(norb: int): + """Test applying fSWAP gate matrix, spinless.""" + + mat01 = np.array([[0, 1], [1, 0]]) + + phase_00 = 1 + phase_11 = -1 + + for i, j in itertools.combinations(range(norb), 2): + for target_orbs in [(i, j), (j, i)]: + assert_has_two_orbital_matrix_spinless( + lambda vec, norb, nelec: ffsim.apply_fswap_gate( + vec, + target_orbs=target_orbs, + norb=norb, + nelec=nelec, + ), + target_orbs=target_orbs, + mat=mat01, + phase_00=phase_00, + phase_11=phase_11, + norb=norb, + )