Hamiltonian
The Hamiltonian class in QDK/Chemistry represents the electronic Hamiltonian operator, which describes the physics of a quantum system.
It contains the one- and two-electron integrals that are essential for quantum chemistry calculations, particularly for active space methods.
Overview
In quantum chemistry, the electronic Hamiltonian is the operator that gives the energy of a system of electrons.
The Hamiltonian class in QDK/Chemistry stores the matrix elements of this operator in the basis of molecular orbitals.
These matrix elements consist of one-electron integrals (representing kinetic energy and electron-nucleus interactions) and two-electron integrals (representing electron-electron repulsion).
Design principles
The Hamiltonian class follows an immutable data model design principle as described in the QDK/Chemistry Design Principles document.
Once properly constructed, the Hamiltonian data is typically not modified during calculations.
This const-correctness approach ensures data integrity throughout computational workflows and prevents accidental modifications of the core quantum system representation.
While setter methods are available for construction and initialization purposes, in normal operation the Hamiltonian object should be treated as immutable after it has been fully populated.
Properties
- One-electron integrals
Matrix of one-electron integrals (h₁)
- Two-electron integrals
Vector of two-electron integrals (h₂) in physicist notation \(\left\langle ij|kl \right\rangle\)
- Core energy
Constant energy term combining nuclear repulsion and inactive orbital contributions
- Inactive Fock matrix
Matrix representing interactions between active and inactive orbitals
- Orbitals
Molecular orbital information for the system (see the Orbitals documentation for detailed information about orbital properties and representations)
- Selected orbital indices
Indices defining the active space orbitals
- Number of electrons
Count of electrons in the active space
Restricted vs. unrestricted Hamiltonians
The Hamiltonian class supports both restricted and unrestricted representations:
Restricted: Uses the same spatial orbitals for alpha and beta electrons. Suitable for closed-shell systems where alpha and beta electrons occupy the same spatial orbitals with opposite spins.
Unrestricted: Allows different spatial orbitals for alpha and beta electrons.
For unrestricted Hamiltonians, the one-electron and two-electron integrals are stored separately for each spin channel:
One-electron integrals: \(h_{\alpha\alpha}\) and \(h_{\beta\beta}\)
Two-electron integrals: \(h_{\alpha\alpha\alpha\alpha}\), \(h_{\alpha\beta\alpha\beta}\), and \(h_{\beta\beta\beta\beta}\)
Usage
The Hamiltonian class is typically used as input to correlation methods such as Configuration Interaction (CI) and Multi-Configuration Self-Consistent Field (MCSCF) calculations.
The HamiltonianConstructor algorithm is the primary tool for generating Hamiltonian objects from molecular data.
Creating a Hamiltonian object
The Hamiltonian object is typically created using the HamiltonianConstructor algorithm (recommended approach for most users), or it can be created directly with the appropriate integral data.
Once properly constructed with all required data, the Hamiltonian object should be considered constant and not modified:
// Load structure from XYZ file
auto structure = Structure::from_xyz_file("../data/water.structure.xyz");
// Run initial SCF to get orbitals
auto scf_solver = ScfSolverFactory::create();
auto [E_scf, wfn] = scf_solver->run(structure, 0, 1);
auto orbitals = wfn->get_orbitals();
// Create a Hamiltonian constructor
// Returns std::shared_ptr<HamiltonianConstructor>
auto hamiltonian_constructor = HamiltonianConstructorFactory::create();
// Construct the Hamiltonian from orbitals
auto hamiltonian = hamiltonian_constructor->run(orbitals);
from pathlib import Path
from qdk_chemistry.algorithms import create
from qdk_chemistry.data import Structure, SpinChannel
# Load a structure from XYZ file
structure = Structure.from_xyz_file(
Path(__file__).parent / "../data/water.structure.xyz"
)
# Run initial SCF to get orbitals
scf_solver = create("scf_solver")
E_scf, wfn = scf_solver.run(
structure, charge=0, spin_multiplicity=1, basis_or_guess="sto-3g"
)
orbitals = wfn.get_orbitals()
# Create a Hamiltonian constructor
hamiltonian_constructor = create("hamiltonian_constructor")
# Construct the Hamiltonian from orbitals
hamiltonian = hamiltonian_constructor.run(orbitals)
Accessing Hamiltonian data
The Hamiltonian class provides methods to access the one- and two-electron integrals and other properties.
In line with its immutable design principle, these methods return const references or copies of the internal data:
Two-electron integral storage and notation
Two-electron integrals in quantum chemistry can be represented using different notations and storage formats. QDK/Chemistry uses the physicist notation by default, but it’s important to understand the different conventions:
- Physicist/Dirac notation \(\left\langle ij|kl \right\rangle\) or \(\left\langle ij|kl \right\rangle\)
Represents the Coulomb interaction where electron 1 occupies orbitals \(i\) and \(k\), while electron 2 occupies orbitals \(j\) and \(l\). This is the default representation in QDK/Chemistry. In this notation, the first index of each pair \((i,k)\) refers to electron 1, and the second index of each pair \((j,l)\) refers to electron 2, following a (1,2,1,2) electron indexing pattern.
- Chemist/Mulliken notation \((ij|kl)\) or \([ij|kl]\)
Represents the Coulomb interaction where electron 1 occupies orbitals \(i\) and \(j\), while electron 2 occupies orbitals \(k\) and \(l\). In this notation, the first pair of indices \((i,j)\) refers to electron 1, and the second pair \((k,l)\) refers to electron 2, following a (1,1,2,2) electron indexing pattern. The symbols differ (parentheses vs square brackets), but the indexing convention is the same.
The relationship between physicist and chemist notation is:
Two-electron integrals with real-valued orbitals possess inherent symmetry properties. From a theoretical perspective, these symmetries can be expressed as:
These permutational symmetries arise from the mathematical properties of the two-electron repulsion integrals.
When accessing specific elements with get_two_body_element(i, j, k, l), the function handles the appropriate index mapping to retrieve the correct value based on the implementation’s storage format.
// Example indices for one- and two-body integral access
size_t i = 0;
size_t j = 1;
size_t k = 2;
size_t l = 3;
// Access one-electron integrals, returns tuple of const Eigen::MatrixXd&
auto [h1_a, h1_b] = hamiltonian.get_one_body_integrals();
// Access two-electron integrals, returns triple of const Eigen::VectorXd&
auto [h2_aaaa, h2_aabb, h2_bbbb] = hamiltonian.get_two_body_integrals();
// Access a specific one-electron integral <ij> (for aa spin channel)
double one_body_element =
hamiltonian.get_one_body_element(i, j, SpinChannel::aa);
// Access a specific two-electron integral <ij|kl> (for aaaa spin channel)
double two_body_element =
hamiltonian.get_two_body_element(i, j, k, l, SpinChannel::aaaa);
// Get core energy (nuclear repulsion + inactive orbital energy), returns double
auto core_energy = hamiltonian.get_core_energy();
// Get orbital data, returns const Orbitals&
const auto& orbitals_access = hamiltonian.get_orbitals();
# Example indices for one- and two-electron integral access
i_int, j_int, k_int, l_int = 0, 1, 2, 3
# Access one-electron integrals (both spin channels)
h1_a, h1_b = hamiltonian.get_one_body_integrals()
# Access two-electron integrals (both spin channels)
h2_aaaa, h2_aabb, h2_bbbb = hamiltonian.get_two_body_integrals()
# Access a specific one-electron integral <ij> (for aa spin channel)
one_body_element = hamiltonian.get_one_body_element(i_int, j_int, SpinChannel.aa)
# Access a specific two-electron integral <ij|kl> (for aaaa spin channel)
two_body_element = hamiltonian.get_two_body_element(
i_int, j_int, k_int, l_int, SpinChannel.aaaa
)
# Get core energy (nuclear repulsion + inactive orbital energy)
core_energy = hamiltonian.get_core_energy()
# Get orbital data
orbitals = hamiltonian.get_orbitals()
Serialization
The Hamiltonian class supports serialization to and from JSON and HDF5 formats.
For detailed information about serialization in QDK/Chemistry, see the Serialization documentation.
Active space Hamiltonian
When constructed with active orbital specifications, the Hamiltonian represents an active space Hamiltonian, which is a projection of the full electronic Hamiltonian into a smaller subspace.
This is essential for tractable multi-configuration calculations.
The HamiltonianConstructor algorithm handles the complex process of generating an appropriate active space Hamiltonian based on your specifications.
Validation methods
The Hamiltonian class provides methods to check the validity and consistency of its data:
// Check if the Hamiltonian data is complete and consistent
// Returns bool
bool valid = hamiltonian.is_valid();
// Check if specific components are available
// All return bool
bool has_one_body = hamiltonian.has_one_body_integrals();
bool has_two_body = hamiltonian.has_two_body_integrals();
bool has_orbitals = hamiltonian.has_orbitals();
Note
This example shows the API pattern. For complete working examples, see the test suite.
# Check if specific components are available
has_one_body = hamiltonian.has_one_body_integrals()
has_two_body = hamiltonian.has_two_body_integrals()
has_orbitals = hamiltonian.has_orbitals()
Further reading
The above examples can be downloaded as complete C++ and Python scripts.
Serialization: Data serialization and deserialization in QDK/Chemistry
Design principles: Design principles for data classes in QDK/Chemistry
Settings: Configuration options for algorithms operating on Hamiltonians