Wavefunction
The Wavefunction class in QDK/Chemistry represents quantum mechanical wavefunctions for molecular systems.
This class provides access to wavefunction coefficients, determinants, reduced density matrices (RDM), orbital entropies [BT15], and other quantum chemical properties.
Overview
A wavefunction in quantum chemistry describes the quantum state of a molecular system.
In QDK/Chemistry, the Wavefunction class encapsulates various wavefunction types, from simple single-determinant Hartree-Fock wavefunctions to complex multi-reference wavefunctions.
The class uses a container-based design where different wavefunction types (Slater determinants, configuration interaction, coupled cluster, etc.) are implemented as specialized container classes, while the main Wavefunction class provides a unified interface.
Mathematical representation
Wavefunctions are represented as linear combinations of determinants:
where \(c_I\) are expansion coefficients and \(|\Phi_I\rangle\) are Slater determinants.
For post-Hartree-Fock methods like coupled cluster, the wavefunction is expressed in terms of cluster operators:
where \(\hat{T} = \hat{T}_1 + \hat{T}_2 + ...\) is the cluster operator and \(|\Phi_0\rangle\) is the reference determinant.
Container types
QDK/Chemistry supports different wavefunction container types for various quantum chemistry methods:
Slater determinant container
Single-determinant wavefunctions (e.g., from Hartree-Fock calculations).
// Use helper function to get orbitals
std::shared_ptr<Orbitals> orbitals = make_minimal_orbitals();
// Create a simple Slater determinant wavefunction for H2 ground state
// 2 electrons in bonding sigma orbital
Configuration det("20");
// Constructor takes single determinant and orbitals as input
auto sd_container =
std::make_unique<SlaterDeterminantContainer>(det, orbitals);
Wavefunction sd_wavefunction(std::move(sd_container));
# Use helper function to get orbitals
orbitals = make_minimal_orbitals()
# Create a simple Slater determinant wavefunction for H2 ground state
# 2 electrons in bonding sigma orbital
det = Configuration("20")
# Constructor takes single determinant and orbitals as input
sd_container = SlaterDeterminantContainer(det, orbitals)
sd_wavefunction = Wavefunction(sd_container)
SCI wavefunction container
Sparse multi-determinant wavefunctions for Selected Configuration Interaction methods.
// Create an SCI wavefunction for H2
// SCI selects only the most important configurations/determinants from the
// full space
std::vector<Configuration> sci_dets = {
Configuration("20"), // Ground state
Configuration("ud"), // Mixed state
Configuration("du") // Mixed state
};
// Coefficients for selected determinants
Eigen::VectorXd sci_coeffs(3);
sci_coeffs << 0.96, 0.15, 0.15;
// Create a SCI wavefunction: requires selected coefficients and determinants,
// as well as orbitals, in constructor
auto sci_container = std::make_unique<SciWavefunctionContainer>(
sci_coeffs, sci_dets, orbitals);
Wavefunction sci_wavefunction(std::move(sci_container));
# Create an SCI wavefunction for H2
# SCI selects only the most important configurations/determinants from the full space
sci_dets = [
Configuration("20"), # both electrons in bonding MO (ground state)
Configuration("du"), # alpha in bonding, beta in antibonding
Configuration("ud"), # beta in bonding, alpha in antibonding
]
# Coefficients for selected determinants
sci_coeffs = np.array([0.96, 0.15, 0.15])
# Create a SCI wavefunction: requires selected coefficients and determinants, as well
# as orbitals, in constructor
sci_container = SciWavefunctionContainer(sci_coeffs, sci_dets, orbitals)
sci_wavefunction = Wavefunction(sci_container)
CAS wavefunction container
A multi-determinant wavefunction from Complete Active Space methods (CASSCF/CASCI).
// Create a CAS wavefunction for H2
// CAS(2,2) = 2 electrons in 2 MOs (bonding and antibonding)
// All possible configurations:
std::vector<Configuration> cas_dets = {
Configuration("20"), // Both electrons in bonding (ground state)
Configuration("ud"), // Alpha in bonding, beta in antibonding
Configuration("du"), // Beta in bonding, alpha in antibonding
Configuration("02") // Both electrons in antibonding
};
// Coefficients
Eigen::VectorXd cas_coeffs(4);
cas_coeffs << 0.95, 0.15, 0.15, 0.05; // Normalized later by the container
// Create a CAS wavefunction : requires all coefficients and determinants, as
// well as orbitals, in constructor
auto cas_container = std::make_unique<CasWavefunctionContainer>(
cas_coeffs, cas_dets, orbitals);
Wavefunction cas_wavefunction(std::move(cas_container));
# Create a CAS wavefunction for H2
# CAS(2,2) = 2 electrons in 2 MOs (bonding and antibonding)
# All possible configurations:
cas_dets = [
Configuration("20"), # both electrons in bonding MO (ground state)
Configuration("ud"), # alpha in bonding, beta in antibonding
Configuration("du"), # beta in bonding, alpha in antibonding
Configuration("02"), # both electrons in antibonding
]
# Coefficients (normalized later by container)
cas_coeffs = np.array([0.95, 0.15, 0.15, 0.05])
# Create a CAS wavefunction: requires all coefficients and determinants,
# as well as orbitals, in constructor
cas_container = CasWavefunctionContainer(cas_coeffs, cas_dets, orbitals)
cas_wavefunction = Wavefunction(cas_container)
MP2 wavefunction container
A wavefunction container for second-order Møller-Plesset perturbation theory (MP2), which stores a reference wavefunction and Hamiltonian. From these, T1 and T2 amplitudes can be computed on demand.
// Create an MP2 wavefunction for H2
// MP2 uses a reference wavefunction and Hamiltonian to compute amplitudes on
// demand
// Use the Slater determinant as reference
auto orbitals_mp2 = make_minimal_orbitals();
auto hamiltonian = make_minimal_hamiltonian(orbitals_mp2);
Configuration ref_det("20");
auto sd_container_mp2 =
std::make_unique<SlaterDeterminantContainer>(ref_det, orbitals_mp2);
auto ref_wavefunction =
std::make_shared<Wavefunction>(std::move(sd_container_mp2));
// Create MP2 container: requires Hamiltonian and reference wavefunction
// Amplitudes are computed lazily when first requested
auto mp2_container =
std::make_unique<MP2Container>(hamiltonian, ref_wavefunction, "mp");
Wavefunction mp2_wavefunction(std::move(mp2_container));
# Create an MP2 wavefunction for H2
# MP2 uses a reference wavefunction and Hamiltonian to compute amplitudes on demand
# Use the Slater determinant as reference
orbitals = make_minimal_orbitals()
hamiltonian = make_minimal_hamiltonian(orbitals)
ref_det = Configuration("20")
sd_container = SlaterDeterminantContainer(ref_det, orbitals)
ref_wavefunction = Wavefunction(sd_container)
# Create MP2 container: requires Hamiltonian and reference wavefunction
# Amplitudes are computed lazily when first requested
mp2_container = MP2Container(hamiltonian, ref_wavefunction, "mp")
mp2_wavefunction = Wavefunction(mp2_container)
Coupled cluster wavefunction container
A coupled cluster wavefunction container that stores T1 and T2 cluster amplitudes along with a reference wavefunction. The container supports reduced density matrices (RDM), which are available if they are provided at construction or computed and stored; otherwise, RDM-related operations are not available.
// Create a coupled cluster wavefunction for H2
// CC uses a reference wavefunction and pre-computed amplitudes
// Use the Slater determinant as reference
auto orbitals_cc = make_minimal_orbitals();
Configuration ref_det_cc("20");
auto sd_container_cc =
std::make_unique<SlaterDeterminantContainer>(ref_det_cc, orbitals_cc);
auto ref_wavefunction_cc =
std::make_shared<Wavefunction>(std::move(sd_container_cc));
// Create example T1 and T2 amplitudes
// T1: occupied-virtual excitations (1 occ × 1 virt = 1 element for H2)
Eigen::VectorXd t1_amplitudes(1);
t1_amplitudes << 0.05;
// T2: occupied-occupied to virtual-virtual excitations
// (1 occ pair × 1 virt pair = 1 element for H2)
Eigen::VectorXd t2_amplitudes(1);
t2_amplitudes << 0.15;
// Create CC container: requires reference wavefunction, orbitals, and
// amplitudes
auto cc_container = std::make_unique<CoupledClusterContainer>(
orbitals_cc, ref_wavefunction_cc, t1_amplitudes, t2_amplitudes);
Wavefunction cc_wavefunction(std::move(cc_container));
# Create a coupled cluster wavefunction for H2
# CC uses a reference wavefunction and pre-computed amplitudes
# Use the Slater determinant as reference
orbitals = make_minimal_orbitals()
ref_det = Configuration("20")
sd_container = SlaterDeterminantContainer(ref_det, orbitals)
ref_wavefunction = Wavefunction(sd_container)
# Create example T1 and T2 amplitudes
# T1: occupied-virtual excitations (1 occ × 1 virt = 1 element for H2)
t1_amplitudes = np.array([0.05])
# T2: occupied-occupied to virtual-virtual excitations
# (1 occ pair × 1 virt pair = 1 element for H2)
t2_amplitudes = np.array([0.15])
# Create CC container: requires reference wavefunction, orbitals, and amplitudes
cc_container = CoupledClusterContainer(
orbitals, ref_wavefunction, t1_amplitudes, t2_amplitudes
)
cc_wavefunction = Wavefunction(cc_container)
Properties
The Wavefunction class provides access to various quantum chemical properties. Availability depends on the specific container type:
Property |
Slater determinant |
Coupled cluster |
|||
|---|---|---|---|---|---|
Coefficients |
✓ |
✓ |
✓ |
✗ |
✗ |
Determinants |
✓ |
✓ |
✓ |
✗ |
✗ |
Electron counts |
✓ |
✓ |
✓ |
✓ |
✓ |
Orbital occupations |
✓ |
✓ |
✓ |
✗ |
✓ |
1-RDMs (spin-dependent) |
✓ |
✓ |
✓ |
✗ |
✓ |
1-RDMs (spin-traced) |
✓ |
✓ |
✓ |
✗ |
✓ |
2-RDMs (spin-dependent) |
✓ |
✓* |
✓* |
✗ |
✓ |
2-RDMs (spin-traced) |
✓ |
✓* |
✓* |
✗ |
✓ |
Orbital entropies |
✓ |
✓* |
✓* |
✗ |
✓ |
T1/T2 amplitudes |
✗ |
✗ |
✗ |
✓ |
✓ |
Overlap calculations |
✗ |
✓ |
✗ |
✗ |
✗ |
Norm calculations |
✓ |
✓ |
✗ |
✗ |
✗ |
Legend:
✓ Available and implemented
✗ Not available (method not implemented)
✓* Implemented and available only if 2-RDMs were provided during construction
Accessing wavefunction data
The Wavefunction class provides methods to access coefficients, determinants, and derived properties:
// Access coefficient(s) and determinant(s) - SD has only one
auto coeffs = sd_wavefunction.get_coefficients();
auto dets = sd_wavefunction.get_active_determinants();
// Get orbital information
auto orbitals_ref = sd_wavefunction.get_orbitals();
// Get electron counts
auto [n_alpha, n_beta] = sd_wavefunction.get_total_num_electrons();
// Get RDMs
auto [rdm1_aa, rdm1_bb] = sd_wavefunction.get_active_one_rdm_spin_dependent();
auto rdm1_total = sd_wavefunction.get_active_one_rdm_spin_traced();
auto [rdm2_aa, rdm2_aabb, rdm2_bbbb] =
sd_wavefunction.get_active_two_rdm_spin_dependent();
auto rdm2_total = sd_wavefunction.get_active_two_rdm_spin_traced();
// Get single orbital entropies
auto entropies = sd_wavefunction.get_single_orbital_entropies();
# Access coefficient(s) and determinant(s) - SD has only one
coeffs = sd_wavefunction.get_coefficients()
dets = sd_wavefunction.get_active_determinants()
# Get orbital information
orbitals_ref = sd_wavefunction.get_orbitals()
# Get electron counts
n_alpha, n_beta = sd_wavefunction.get_total_num_electrons()
# Get RDMs
rdm1_aa, rdm1_bb = sd_wavefunction.get_active_one_rdm_spin_dependent()
rdm1_total = sd_wavefunction.get_active_one_rdm_spin_traced()
rdm2_aaaa, rdm2_aabb, rdm2_bbbb = sd_wavefunction.get_active_two_rdm_spin_dependent()
rdm2_total = sd_wavefunction.get_active_two_rdm_spin_traced()
# Get single orbital entropies
entropies = sd_wavefunction.get_single_orbital_entropies()
Accessing cluster amplitudes
For MP2 and coupled cluster wavefunctions, one can access T1 and T2 cluster amplitudes.
// Access T1 and T2 amplitudes from MP2 and CC containers
// MP2
// Get the container back from wfn
const auto& mp2_container_ref =
mp2_wavefunction.get_container<MP2Container>();
// Amplitudes are lazily evaluated on first call then cached
auto [t2_abab_mp2, t2_aaaa_mp2, t2_bbbb_mp2] =
mp2_container_ref->get_t2_amplitudes();
// CC
const auto& cc_container_ref =
cc_wavefunction.get_container<CoupledClusterContainer>();
// Amplitudes are stored already from construction
auto [t1_aa, t1_bb] = cc_container_ref->get_t1_amplitudes();
auto [t2_abab_cc, t2_aaaa_cc, t2_bbbb_cc] =
cc_container_ref->get_t2_amplitudes();
# Access T1 and T2 amplitudes from MP2 and CC containers
# For MP2 - amplitudes computed on demand
t2_abab_mp2, t2_aaaa_mp2, t2_bbbb_mp2 = (
mp2_wavefunction.get_container().get_t2_amplitudes()
)
# For CC - amplitudes stored during construction
t1_aa, t1_bb = cc_wavefunction.get_container().get_t1_amplitudes()
t2_abab_cc, t2_aaaa_cc, t2_bbbb_cc = cc_wavefunction.get_container().get_t2_amplitudes()
Further reading
Settings: Configuration settings for algorithms
Serialization: Data serialization and deserialization
Active space methods: Active space selection from wavefunctions