Hamiltonians

The hamiltonians subpackage provides abstract and concrete Hamiltonian operators for molecular and spin systems.

Base Classes

class qvartools.hamiltonians.hamiltonian.Hamiltonian(num_sites, local_dim=2)[source]

Bases: ABC

Abstract base class for all Hamiltonian operators.

Every subclass must implement diagonal_element() and get_connections(). The base class provides dense/sparse matrix construction, exact diagonalisation, and configuration bookkeeping.

Parameters:
  • num_sites (int) – Number of lattice / orbital sites.

  • local_dim (int, optional) – Dimension of the local Hilbert space on each site (default 2 for spin-½ / qubit).

num_sites

Number of sites.

Type:

int

local_dim

Local Hilbert-space dimension.

Type:

int

hilbert_dim

Full Hilbert-space dimension local_dim ** num_sites.

Type:

int

Notes

This is an abstract class and cannot be instantiated directly. Subclasses must implement diagonal_element() and get_connections().

Examples

See HeisenbergHamiltonian or TransverseFieldIsing for concrete implementations.

Methods

diagonal_element(config)

Return the diagonal matrix element ⟨config|H|config⟩.

exact_ground_state([device])

Compute the exact ground state via full diagonalisation.

get_connections(config)

Return all states connected to config by the Hamiltonian.

ground_state_sparse([k, device])

Compute the lowest k eigenstates via sparse diagonalisation.

matrix_element(config_i, config_j)

Compute a single matrix element ⟨config_i|H|config_j⟩.

matrix_elements(configs_bra, configs_ket)

Compute the matrix H_{ij} = ⟨bra_i|H|ket_j⟩.

to_dense([device])

Build the full dense Hamiltonian matrix.

to_sparse([device])

Build a SciPy CSR sparse Hamiltonian matrix.

abstractmethod diagonal_element(config)[source]

Return the diagonal matrix element ⟨config|H|config⟩.

Parameters:

config (torch.Tensor) – Basis-state configuration, shape (num_sites,) with integer entries in [0, local_dim).

Returns:

Scalar tensor (shape ()) with the diagonal element.

Return type:

torch.Tensor

abstractmethod get_connections(config)[source]

Return all states connected to config by the Hamiltonian.

Parameters:

config (torch.Tensor) – Basis-state configuration, shape (num_sites,).

Return type:

Tuple[Tensor, Tensor]

Returns:

  • connected_configs (torch.Tensor) – Off-diagonal connected configurations, shape (n_conn, num_sites).

  • matrix_elements (torch.Tensor) – Corresponding matrix elements ⟨connected|H|config⟩, shape (n_conn,).

matrix_element(config_i, config_j)[source]

Compute a single matrix element ⟨config_i|H|config_j⟩.

Parameters:
  • config_i (torch.Tensor) – Bra configuration, shape (num_sites,).

  • config_j (torch.Tensor) – Ket configuration, shape (num_sites,).

Returns:

Scalar tensor with the matrix element.

Return type:

torch.Tensor

Examples

>>> from qvartools.hamiltonians import TransverseFieldIsing
>>> H = TransverseFieldIsing(num_spins=2, V=1.0, h=0.5)
>>> import torch
>>> H.matrix_element(torch.tensor([0, 0]), torch.tensor([0, 0]))
tensor(...)
matrix_elements(configs_bra, configs_ket)[source]

Compute the matrix H_{ij} = ⟨bra_i|H|ket_j⟩.

Parameters:
  • configs_bra (torch.Tensor) – Bra configurations, shape (M, num_sites).

  • configs_ket (torch.Tensor) – Ket configurations, shape (N, num_sites).

Returns:

Matrix of shape (M, N) with dtype float64.

Return type:

torch.Tensor

to_dense(device='cpu')[source]

Build the full dense Hamiltonian matrix.

Parameters:

device (str, optional) – Torch device for the output tensor (default "cpu").

Returns:

Dense matrix of shape (hilbert_dim, hilbert_dim), dtype float64.

Return type:

torch.Tensor

Warns:

UserWarning – If num_sites > 16, since the matrix may be very large.

Examples

>>> from qvartools.hamiltonians import TransverseFieldIsing
>>> H = TransverseFieldIsing(num_spins=2, V=1.0, h=0.5)
>>> H.to_dense().shape
torch.Size([4, 4])
to_sparse(device='cpu')[source]

Build a SciPy CSR sparse Hamiltonian matrix.

Parameters:

device (str, optional) – Torch device used during construction (default "cpu"). The returned sparse matrix always lives on the CPU (NumPy).

Returns:

Sparse matrix of shape (hilbert_dim, hilbert_dim).

Return type:

scipy.sparse.csr_matrix

exact_ground_state(device='cpu')[source]

Compute the exact ground state via full diagonalisation.

Parameters:

device (str, optional) – Torch device (default "cpu").

Return type:

Tuple[float, Tensor]

Returns:

  • energy (float) – Ground-state energy (lowest eigenvalue).

  • state (torch.Tensor) – Ground-state eigenvector, shape (hilbert_dim,).

Examples

>>> from qvartools.hamiltonians import TransverseFieldIsing
>>> H = TransverseFieldIsing(num_spins=2, V=1.0, h=0.5)
>>> energy, state = H.exact_ground_state()
>>> isinstance(energy, float)
True
ground_state_sparse(k=1, device='cpu')[source]

Compute the lowest k eigenstates via sparse diagonalisation.

Uses scipy.sparse.linalg.eigsh() with which='SA'.

Parameters:
  • k (int, optional) – Number of lowest eigenstates to compute (default 1).

  • device (str, optional) – Torch device used during sparse matrix construction (default "cpu").

Return type:

Tuple[ndarray, ndarray]

Returns:

  • eigenvalues (np.ndarray) – Lowest k eigenvalues, shape (k,).

  • eigenvectors (np.ndarray) – Corresponding eigenvectors, shape (hilbert_dim, k).

class qvartools.hamiltonians.pauli_string.PauliString(paulis, coefficient=1.0)[source]

Bases: object

A single Pauli string (tensor product of I, X, Y, Z) with coefficient.

Parameters:
  • paulis (list of str) – One-character Pauli labels for each qubit. Allowed values are "I", "X", "Y", "Z".

  • coefficient (complex, optional) – Multiplicative coefficient (default 1.0).

paulis

The Pauli labels.

Type:

list of str

coefficient

The coefficient.

Type:

complex

num_qubits

Number of qubits the string acts on.

Type:

int

Raises:

ValueError – If any label in paulis is not in {"I", "X", "Y", "Z"}.

Parameters:

Examples

>>> ps = PauliString(["X", "Z", "I"], coefficient=0.5)
>>> new_config, phase = ps.apply(torch.tensor([0, 1, 0]))

Methods

apply(config)

Apply this Pauli string to a computational-basis state.

is_diagonal()

Return True if the Pauli string is diagonal (only I and Z).

apply(config)[source]

Apply this Pauli string to a computational-basis state.

The computational basis is |0⟩, |1⟩ per qubit. The rules are:

  • I — identity, no change.

  • X — bit flip: |0⟩ |1⟩, |1⟩ |0⟩.

  • Y — bit flip with phase: |0⟩ i|1⟩, |1⟩ −i|0⟩.

  • Z — phase flip: |0⟩ |0⟩, |1⟩ −|1⟩.

Parameters:

config (torch.Tensor) – Input configuration, shape (num_qubits,) with entries in {0, 1}.

Return type:

Tuple[Optional[Tensor], complex]

Returns:

  • new_config (torch.Tensor or None) – The resulting configuration after applying the string. None is never returned for valid inputs (kept in the signature for forward-compatibility with annihilation operators).

  • phase (complex) – The accumulated phase (including self.coefficient).

Examples

>>> ps = PauliString(["X", "Z"], coefficient=1.0)
>>> new_cfg, phase = ps.apply(torch.tensor([0, 1]))
>>> new_cfg
tensor([1, 1])
>>> phase
(-1+0j)
is_diagonal()[source]

Return True if the Pauli string is diagonal (only I and Z).

Return type:

bool

Examples

>>> PauliString(["Z", "I", "Z"]).is_diagonal()
True
>>> PauliString(["X", "I"]).is_diagonal()
False

Molecular Hamiltonians

class qvartools.hamiltonians.molecular.hamiltonian.MolecularHamiltonian(integrals, device='cpu')[source]

Bases: Hamiltonian

Second-quantised molecular Hamiltonian with Jordan–Wigner mapping.

Spin-orbitals are ordered as [alpha_0, alpha_1, ..., alpha_{n-1}, beta_0, beta_1, ..., beta_{n-1}]. A configuration is a binary occupation vector of length 2 * n_orbitals.

Parameters:
  • integrals (MolecularIntegrals) – Molecular integrals (one-body, two-body, metadata).

  • device (str, optional) – Torch device for tensor storage (default "cpu").

integrals

The stored integrals.

Type:

MolecularIntegrals

n_orb

Number of spatial orbitals.

Type:

int

E_nuc

Nuclear repulsion energy.

Type:

float

device

Torch device string.

Type:

str

Examples

>>> import numpy as np
>>> from qvartools.hamiltonians import MolecularIntegrals, MolecularHamiltonian
>>> h1 = np.diag([-1.0, -0.5])
>>> h2 = np.zeros((2, 2, 2, 2))
>>> mi = MolecularIntegrals(h1, h2, 0.7, 2, 2, 1, 1)
>>> ham = MolecularHamiltonian(mi)
>>> ham.num_sites
4
Attributes:
h1e

One-electron integral tensor.

h2e

Two-electron integral tensor.

n_alpha

Number of alpha electrons.

n_beta

Number of beta electrons.

n_orbitals

Number of spatial orbitals.

Parameters:

Methods

diagonal_element(config)

Compute the diagonal element ⟨config|H|config⟩.

diagonal_elements_batch(configs)

Compute diagonal matrix elements for a batch of configurations.

exact_ground_state([device])

Compute the exact ground state via full diagonalisation.

fci_energy()

Compute the Full-CI energy using PySCF's native FCI solver.

get_connections(config)

Return all off-diagonal connected states via Slater--Condon rules.

get_connections_vectorized_batch(configs)

Compute connections for a batch of configs (GPU-vectorised).

get_hf_state()

Return the Hartree--Fock ground-state configuration.

ground_state_sparse([k, device])

Compute the lowest k eigenstates via sparse diagonalisation.

matrix_element(config_i, config_j)

Compute a single matrix element ⟨config_i|H|config_j⟩.

matrix_elements(configs_bra, configs_ket)

Compute H_{ij} = ⟨bra_i|H|ket_j⟩ using hash-based lookup.

matrix_elements_fast(configs)

Build a Hermitian projected Hamiltonian matrix efficiently.

to_dense([device])

Build the full dense Hamiltonian matrix.

to_sparse([device])

Build a SciPy CSR sparse Hamiltonian matrix.

diagonal_elements_batch(configs)[source]

Compute diagonal matrix elements for a batch of configurations.

Uses vectorised einsum operations for the one-body, Coulomb, and exchange contributions.

Parameters:

configs (torch.Tensor) – Occupation vectors, shape (batch, num_sites), dtype int64 or float64.

Returns:

Diagonal elements, shape (batch,), dtype float64.

Return type:

torch.Tensor

diagonal_element(config)[source]

Compute the diagonal element ⟨config|H|config⟩.

Parameters:

config (torch.Tensor) – Occupation vector, shape (num_sites,).

Returns:

Scalar tensor with the diagonal energy.

Return type:

torch.Tensor

get_connections(config)[source]

Return all off-diagonal connected states via Slater–Condon rules.

Computes single and double excitations with proper Jordan–Wigner signs. Uses Numba kernels when available, otherwise falls back to a pure-Python implementation.

Parameters:

config (torch.Tensor) – Occupation vector, shape (num_sites,).

Return type:

Tuple[Tensor, Tensor]

Returns:

  • connected_configs (torch.Tensor) – Shape (n_conn, num_sites), dtype int64.

  • matrix_elements (torch.Tensor) – Shape (n_conn,), dtype float64.

get_connections_vectorized_batch(configs)[source]

Compute connections for a batch of configs (GPU-vectorised).

This method processes each configuration independently but uses GPU tensors throughout to avoid CPU↔GPU transfers within each configuration’s excitation enumeration.

Parameters:

configs (torch.Tensor) – Occupation vectors, shape (batch, num_sites).

Return type:

Tuple[List[Tensor], List[Tensor]]

Returns:

  • all_connected (list of torch.Tensor) – One tensor per batch element, each of shape (n_conn_i, num_sites).

  • all_elements (list of torch.Tensor) – One tensor per batch element, each of shape (n_conn_i,).

matrix_elements(configs_bra, configs_ket)[source]

Compute H_{ij} = ⟨bra_i|H|ket_j⟩ using hash-based lookup.

This override of the base-class method uses integer hashing to avoid the O(n_conn * M) inner loop for matching bra states. Uses vectorised batch hashing for both bra and connected configs.

Parameters:
  • configs_bra (torch.Tensor) – Bra configurations, shape (M, num_sites).

  • configs_ket (torch.Tensor) – Ket configurations, shape (N, num_sites).

Returns:

Matrix of shape (M, N), dtype float64.

Return type:

torch.Tensor

matrix_elements_fast(configs)[source]

Build a Hermitian projected Hamiltonian matrix efficiently.

Builds only the lower triangle via hash-based matching then mirrors to the upper triangle, guaranteeing exact symmetry.

Uses fully vectorised hash lookups: torch.searchsorted on sorted basis hashes replaces the per-connection Python loop, giving O(n_conn * log(n_configs)) matching per ket instead of O(n_conn) Python dict lookups.

Parameters:

configs (torch.Tensor) – Basis configurations, shape (n_configs, num_sites).

Returns:

Hermitian matrix of shape (n_configs, n_configs), dtype float64.

Return type:

torch.Tensor

Raises:

MemoryError – If n_configs > 10000 (dense matrix would be too large).

property n_orbitals: int

Number of spatial orbitals.

property n_alpha: int

Number of alpha electrons.

property n_beta: int

Number of beta electrons.

property h1e: Tensor

One-electron integral tensor.

property h2e: Tensor

Two-electron integral tensor.

get_hf_state()[source]

Return the Hartree–Fock ground-state configuration.

Alpha electrons fill the lowest alpha spin-orbitals, beta electrons fill the lowest beta spin-orbitals.

Returns:

Binary occupation vector, shape (num_sites,), dtype int64.

Return type:

torch.Tensor

Examples

>>> import numpy as np
>>> from qvartools.hamiltonians import MolecularIntegrals, MolecularHamiltonian
>>> mi = MolecularIntegrals(np.eye(2), np.zeros((2,2,2,2)), 0.0, 2, 2, 1, 1)
>>> ham = MolecularHamiltonian(mi)
>>> ham.get_hf_state()
tensor([1, 0, 0, 1])
fci_energy()[source]

Compute the Full-CI energy using PySCF’s native FCI solver.

Uses PySCF’s Davidson-iteration FCI solver in the compressed alpha/beta string representation, which is vastly faster than building the full 2^n dense matrix. Falls back to the GPU FCI solver if CuPy is available, or to brute-force dense diagonalisation for very small systems (≤ 8 qubits).

Returns:

The ground-state energy in Hartree.

Return type:

float

Raises:

RuntimeError – If PySCF is not installed and num_sites > 20 (dense diagonalisation would be intractable).

class qvartools.hamiltonians.integrals.MolecularIntegrals(h1e, h2e, nuclear_repulsion, n_electrons, n_orbitals, n_alpha, n_beta)[source]

Bases: object

Container for molecular one- and two-electron integrals.

All arrays use the spatial-orbital indexing convention produced by PySCF’s ao2mo module.

Parameters:
  • h1e (np.ndarray) – One-electron integrals, shape (n_orb, n_orb), dtype float64.

  • h2e (np.ndarray) – Two-electron integrals in chemist’s notation (pq|rs), shape (n_orb, n_orb, n_orb, n_orb), dtype float64.

  • nuclear_repulsion (float) – Nuclear repulsion energy in Hartree.

  • n_electrons (int) – Total number of electrons.

  • n_orbitals (int) – Number of spatial orbitals.

  • n_alpha (int) – Number of alpha (spin-up) electrons.

  • n_beta (int) – Number of beta (spin-down) electrons.

Raises:

ValueError – If array shapes are inconsistent with n_orbitals or if n_alpha + n_beta != n_electrons.

Examples

>>> import numpy as np
>>> h1 = np.zeros((2, 2))
>>> h2 = np.zeros((2, 2, 2, 2))
>>> mi = MolecularIntegrals(h1, h2, 0.7, 2, 2, 1, 1)
>>> mi.n_orbitals
2
h1e: ndarray
h2e: ndarray
nuclear_repulsion: float
n_electrons: int
n_orbitals: int
n_alpha: int
n_beta: int
qvartools.hamiltonians.integrals.compute_molecular_integrals(geometry, basis='sto-3g', charge=0, spin=0)[source]

Run RHF with PySCF and extract molecular integrals.

Parameters:
  • geometry (list of (str, (float, float, float))) – Molecular geometry. Each element is (atom_symbol, (x, y, z)) with coordinates in Angstroms.

  • basis (str, optional) – Gaussian basis set name (default "sto-3g").

  • charge (int, optional) – Net charge of the molecule (default 0).

  • spin (int, optional) – Spin multiplicity minus one, i.e. 2S (default 0 for singlet).

Returns:

Integrals and metadata needed by MolecularHamiltonian.

Return type:

MolecularIntegrals

Raises:

Examples

>>> geometry = [("H", (0.0, 0.0, 0.0)), ("H", (0.0, 0.0, 0.74))]
>>> mi = compute_molecular_integrals(geometry, basis="sto-3g")
>>> mi.n_orbitals
2

Spin Hamiltonians

class qvartools.hamiltonians.spin.heisenberg.HeisenbergHamiltonian(num_spins, Jx=1.0, Jy=1.0, Jz=1.0, h_x=0.0, h_y=0.0, h_z=0.0, periodic=True)[source]

Bases: Hamiltonian

Anisotropic Heisenberg (XYZ) model on a 1-D chain.

\[H = \sum_{\langle i,j \rangle} \bigl( J_x\, S^x_i S^x_j + J_y\, S^y_i S^y_j + J_z\, S^z_i S^z_j \bigr) + \sum_i \bigl( h^x_i S^x_i + h^y_i S^y_i + h^z_i S^z_i \bigr)\]

where S operators are spin-½ operators (with eigenvalues ±½) and the sum runs over nearest-neighbour pairs on a chain.

Convention: config[i] = 0 → spin-up (S^z = +½), config[i] = 1 → spin-down (S^z = −½).

Parameters:
  • num_spins (int) – Number of spin-½ sites.

  • Jx (float, optional) – XX coupling constant (default 1.0).

  • Jy (float, optional) – YY coupling constant (default 1.0).

  • Jz (float, optional) – ZZ coupling constant (default 1.0).

  • h_x (float or array-like, optional) – External field in the x-direction (default 0).

  • h_y (float or array-like, optional) – External field in the y-direction (default 0).

  • h_z (float or array-like, optional) – External field in the z-direction (default 0).

  • periodic (bool, optional) – Whether to use periodic boundary conditions (default True).

Examples

>>> import torch
>>> H = HeisenbergHamiltonian(num_spins=4, Jx=1.0, Jy=1.0, Jz=1.0)
>>> config = torch.tensor([0, 1, 0, 1])
>>> H.diagonal_element(config)
tensor(...)
>>> energy, state = H.exact_ground_state()

Methods

diagonal_element(config)

Compute ⟨config|H|config⟩ (ZZ interaction + Z field).

diagonal_elements_batch(configs)

Vectorised diagonal elements for a batch of configurations.

exact_ground_state([device])

Compute the exact ground state via full diagonalisation.

get_connections(config)

Return off-diagonal connected states (XX + YY exchange, X/Y fields).

ground_state_sparse([k, device])

Compute the lowest k eigenstates via sparse diagonalisation.

matrix_element(config_i, config_j)

Compute a single matrix element ⟨config_i|H|config_j⟩.

matrix_elements(configs_bra, configs_ket)

Compute the matrix H_{ij} = ⟨bra_i|H|ket_j⟩.

to_dense([device])

Build the full dense Hamiltonian matrix.

to_sparse([device])

Build a SciPy CSR sparse Hamiltonian matrix.

diagonal_element(config)[source]

Compute ⟨config|H|config⟩ (ZZ interaction + Z field).

Parameters:

config (torch.Tensor) – Spin configuration, shape (num_sites,), entries in {0, 1}.

Returns:

Scalar diagonal energy.

Return type:

torch.Tensor

diagonal_elements_batch(configs)[source]

Vectorised diagonal elements for a batch of configurations.

Parameters:

configs (torch.Tensor) – Shape (batch, num_sites).

Returns:

Shape (batch,), dtype float64.

Return type:

torch.Tensor

get_connections(config)[source]

Return off-diagonal connected states (XX + YY exchange, X/Y fields).

The XX + YY interaction flips a pair of anti-aligned neighbours: |↑↓⟩ |↓↑⟩ with amplitude 0.5 * (Jx + Jy) / 2. When Jx Jy the |↑↓⟩ |↓↑⟩ and |↓↑⟩ |↑↓⟩ amplitudes differ.

The transverse fields h_x and h_y produce single-site flips.

Parameters:

config (torch.Tensor) – Shape (num_sites,).

Return type:

Tuple[Tensor, Tensor]

Returns:

class qvartools.hamiltonians.spin.tfim.TransverseFieldIsing(num_spins, V=1.0, h=1.0, L=1, periodic=True)[source]

Bases: Hamiltonian

Transverse-field Ising model with tuneable interaction range.

\[H = -V \sum_{i} \sum_{l=1}^{L} S^z_i\, S^z_{i+l} - h \sum_i S^x_i\]

where L controls the range of the ZZ interaction and the sum over i respects boundary conditions.

Convention: config[i] = 0 → spin-up (S^z = +½), config[i] = 1 → spin-down (S^z = −½).

Parameters:
  • num_spins (int) – Number of spin sites.

  • V (float, optional) – ZZ interaction strength (default 1.0).

  • h (float, optional) – Transverse field strength (default 1.0).

  • L (int, optional) – Range of the ZZ interaction in lattice units (default 1 for nearest-neighbour).

  • periodic (bool, optional) – Periodic boundary conditions (default True).

Examples

>>> import torch
>>> H = TransverseFieldIsing(num_spins=4, V=1.0, h=0.5)
>>> config = torch.tensor([0, 1, 0, 1])
>>> H.diagonal_element(config)
tensor(...)
>>> energy, state = H.exact_ground_state()

Methods

diagonal_element(config)

Compute ⟨config|H|config⟩ (ZZ interaction with range L).

exact_ground_state([device])

Compute the exact ground state via full diagonalisation.

get_connections(config)

Return off-diagonal connected states from the transverse field.

ground_state_sparse([k, device])

Compute the lowest k eigenstates via sparse diagonalisation.

matrix_element(config_i, config_j)

Compute a single matrix element ⟨config_i|H|config_j⟩.

matrix_elements(configs_bra, configs_ket)

Compute the matrix H_{ij} = ⟨bra_i|H|ket_j⟩.

to_dense([device])

Build the full dense Hamiltonian matrix.

to_sparse([device])

Build a SciPy CSR sparse Hamiltonian matrix.

diagonal_element(config)[source]

Compute ⟨config|H|config⟩ (ZZ interaction with range L).

Parameters:

config (torch.Tensor) – Shape (num_sites,).

Returns:

Scalar diagonal energy.

Return type:

torch.Tensor

get_connections(config)[source]

Return off-diagonal connected states from the transverse field.

The transverse field -h S^x_i flips each spin individually with amplitude -h/2.

Parameters:

config (torch.Tensor) – Shape (num_sites,).

Return type:

Tuple[Tensor, Tensor]

Returns:

Jordan-Wigner Mapping

jordan_wigner — Jordan–Wigner sign computation kernels

Provides Numba-accelerated (with pure-Python fallback) functions for computing the fermionic sign factors arising from the Jordan–Wigner transformation for single and double excitations.

The module also exports the _HAS_NUMBA flag and the njit shim so that downstream modules can decorate their own kernels consistently.

qvartools.hamiltonians.molecular.jordan_wigner.numba_jw_sign_single(config, p, q)

Compute the Jordan–Wigner sign for a single excitation a†_p a_q.

The sign is (-1)^(number of occupied orbitals strictly between p and q).

Parameters:
  • config (np.ndarray) – Occupation vector, shape (num_sites,), dtype int.

  • p (int) – Creation orbital index.

  • q (int) – Annihilation orbital index.

Returns:

+1 or -1.

Return type:

int

Notes

The Jordan–Wigner sign arises from anti-commuting the fermionic operator past all occupied orbitals between positions p and q:

\[ext{sign} = (-1)^{\sum_{k=\min(p,q)+1}^{\max(p,q)-1} n_k}\]

Examples

>>> import numpy as np
>>> config = np.array([1, 1, 0, 1])
>>> numba_jw_sign_single(config, 0, 3)
-1
qvartools.hamiltonians.molecular.jordan_wigner.numba_jw_sign_double(config, p, r, q, s)

Compute the Jordan–Wigner sign for a double excitation a†_p a†_r a_s a_q.

The operator ordering is right-to-left: first annihilate q, then s, then create r, then create p. Each step accumulates a sign from the occupied orbitals it must anti-commute past (on the current state of the configuration after previous operations).

Parameters:
  • config (np.ndarray) – Occupation vector, shape (num_sites,), dtype int.

  • p (int) – First creation orbital.

  • r (int) – Second creation orbital.

  • q (int) – First annihilation orbital.

  • s (int) – Second annihilation orbital.

Returns:

+1 or -1.

Return type:

int

Notes

Each fermionic operator is translated into a Jordan–Wigner string that requires a Z-chain parity count. The four operations are applied sequentially on a mutable copy of the configuration to capture the parity of occupied orbitals at each step.

Examples

>>> import numpy as np
>>> config = np.array([1, 1, 0, 0])
>>> numba_jw_sign_double(config, 2, 3, 0, 1)
1

Slater-Condon Rules

slater_condon — Slater–Condon excitation kernels

Provides Numba-accelerated (with pure-Python fallback) functions for enumerating single and double excitations from a given occupation configuration using Slater–Condon rules together with Jordan–Wigner sign factors.

Functions

_numba_single_excitations

All single excitations from a configuration.

_numba_double_excitations

All double excitations from a configuration.

numba_get_connections

Combined single + double excitations.

qvartools.hamiltonians.molecular.slater_condon.numba_get_connections(config, n_orb, J_single, K_single, h1e, h2e, num_sites)

Compute all single and double excitations via Numba.

Parameters:
  • config (np.ndarray) – Occupation vector, shape (num_sites,).

  • n_orb (int) – Number of spatial orbitals.

  • J_single (np.ndarray) – Coulomb integrals for singles, shape (n_orb, n_orb, n_orb).

  • K_single (np.ndarray) – Exchange integrals for singles, shape (n_orb, n_orb, n_orb).

  • h1e (np.ndarray) – One-electron integrals, shape (n_orb, n_orb).

  • h2e (np.ndarray) – Two-electron integrals, shape (n_orb, n_orb, n_orb, n_orb).

  • num_sites (int) – Total number of spin-orbitals.

Return type:

Tuple[ndarray, ndarray]

Returns:

  • configs (np.ndarray) – Connected configurations, shape (n_conn, num_sites).

  • elements (np.ndarray) – Matrix elements, shape (n_conn,).