Namespace qdk::chemistry::utils::model_hamiltonians
-
namespace model_hamiltonians
Functions
-
template<typename EpsT, typename TT, typename UT>
inline qdk::chemistry::data::Hamiltonian create_hubbard_hamiltonian(const qdk::chemistry::data::LatticeGraph &lattice, EpsT &&epsilon_in, TT &&t_in, UT &&U_in) Create a Hubbard model Hamiltonian.
- Template Parameters:
EpsT – double or Eigen::VectorXd
TT – double or Eigen::MatrixXd
UT – double or Eigen::VectorXd
- Parameters:
lattice – Symmetric lattice graph defining the connectivity.
epsilon_in – On-site orbital energies. Scalar or VectorXd of size n.
t_in – Hopping integrals. Scalar or n x n MatrixXd.
U_in – On-site Coulomb repulsion. Scalar or VectorXd of size n.
- Returns:
Hamiltonian for the Hubbard model.
-
template<typename EpsT, typename TT>
inline qdk::chemistry::data::Hamiltonian create_huckel_hamiltonian(const qdk::chemistry::data::LatticeGraph &lattice, EpsT &&epsilon_in, TT &&t_in) Create a Hückel model Hamiltonian.
- Template Parameters:
EpsT – double or Eigen::VectorXd
TT – double or Eigen::MatrixXd
- Parameters:
lattice – Symmetric lattice graph defining the connectivity.
epsilon_in – On-site orbital energies. Scalar or VectorXd of size n.
t_in – Hopping integrals. Scalar or n x n MatrixXd.
- Returns:
Hamiltonian for the Hückel model.
-
template<typename EpsT, typename TT, typename UT, typename VT, typename ZT>
inline qdk::chemistry::data::Hamiltonian create_ppp_hamiltonian(const qdk::chemistry::data::LatticeGraph &lattice, EpsT &&epsilon_in, TT &&t_in, UT &&U_in, VT &&V_in, ZT &&z_in) Create a Pariser-Parr-Pople (PPP) model Hamiltonian.
- Template Parameters:
EpsT – double or Eigen::VectorXd
TT – double or Eigen::MatrixXd
UT – double or Eigen::VectorXd
VT – double or Eigen::MatrixXd
ZT – double or Eigen::VectorXd
- Parameters:
lattice – Symmetric lattice graph defining the connectivity.
epsilon_in – On-site orbital energies. Scalar or VectorXd of size n.
t_in – Hopping integrals. Scalar or n x n MatrixXd.
U_in – On-site Coulomb repulsion. Scalar or VectorXd of size n.
V_in – Intersite Coulomb interaction matrix. Scalar or n x n MatrixXd.
z_in – Effective core charges. Scalar or VectorXd of size n.
- Returns:
Hamiltonian for the PPP model.
-
template<typename UT, typename RT>
inline Eigen::MatrixXd mataga_nishimoto_potential(const qdk::chemistry::data::LatticeGraph &lattice, UT &&U, RT &&R, double epsilon_r = 1.0, bool nearest_neighbor_only = false) Compute the Mataga-Nishimoto intersite potential matrix.
V_ij = U_ij / (1 + U_ij * epsilon_r * R_ij)
where U_ij = sqrt(U_i * U_j) is the geometric mean of on-site parameters.
All parameters should be in atomic units (Hartree for U, Bohr for R).
- Template Parameters:
UT – double or Eigen::VectorXd
RT – double or Eigen::MatrixXd
- Parameters:
lattice – Lattice graph (used for the number of sites).
U – On-site Coulomb parameter(s) in Hartree. Scalar or VectorXd.
R – Intersite distances in Bohr. Scalar or n x n MatrixXd.
epsilon_r – Relative permittivity (dimensionless, default 1.0).
nearest_neighbor_only – If true, restrict to lattice-connected pairs (default false).
- Returns:
n x n symmetric MatrixXd of Mataga-Nishimoto potential values in Hartree.
-
template<typename UT, typename RT>
inline Eigen::MatrixXd ohno_potential(const qdk::chemistry::data::LatticeGraph &lattice, UT &&U, RT &&R, double epsilon_r = 1.0, bool nearest_neighbor_only = false) Compute the Ohno intersite potential matrix.
V_ij = U_ij / sqrt(1 + (U_ij * epsilon_r * R_ij)^2)
where U_ij = sqrt(U_i * U_j) is the geometric mean of on-site parameters.
All parameters should be in atomic units (Hartree for U, Bohr for R).
- Template Parameters:
UT – double or Eigen::VectorXd
RT – double or Eigen::MatrixXd
- Parameters:
lattice – Lattice graph (used for the number of sites).
U – On-site Coulomb parameter(s) in Hartree. Scalar or VectorXd.
R – Intersite distances in Bohr. Scalar or n x n MatrixXd.
epsilon_r – Relative permittivity (dimensionless, default 1.0).
nearest_neighbor_only – If true, restrict to lattice-connected pairs (default false).
- Returns:
n x n symmetric MatrixXd of Ohno potential values in Hartree.
-
template<typename UT, typename RT, typename PotentialFunc>
inline Eigen::MatrixXd pairwise_potential(const qdk::chemistry::data::LatticeGraph &lattice, UT &&U_in, RT &&R_in, PotentialFunc &&func, bool nearest_neighbor_only = false) Compute a symmetric pairwise potential matrix from a custom formula.
For each unique pair (i < j), computes the geometric mean U_ij = sqrt(U_i * U_j), reads R_ij, and evaluates func(i, j, U_ij, R_ij). The result is stored symmetrically: V(i,j) = V(j,i).
When nearest_neighbor_only is true, only pairs connected by a lattice edge are evaluated; all other entries remain zero.
- Template Parameters:
UT – double or Eigen::VectorXd — on-site Coulomb parameter(s).
RT – double or Eigen::MatrixXd — intersite distances.
PotentialFunc – Callable with signature (int i, int j, double Uij, double Rij) -> double.
- Parameters:
lattice – Lattice graph defining the connectivity and number of sites.
U_in – On-site Coulomb parameter(s). Scalar or VectorXd of size n.
R_in – Distance matrix. Scalar or n x n MatrixXd.
func – Potential formula to evaluate for each pair.
nearest_neighbor_only – If true, restrict to lattice-connected pairs (default false).
- Throws:
std::invalid_argument – if U size or R dimensions mismatch.
- Returns:
n x n symmetric MatrixXd of pairwise potential values.
-
namespace detail
Functions
-
template<typename EpsT, typename TT, typename UT>
inline std::tuple<Eigen::SparseMatrix<double>, qdk::chemistry::data::SparseHamiltonianContainer::TwoBodyMap> _build_hubbard_integrals(const qdk::chemistry::data::LatticeGraph &lattice, EpsT &&epsilon_in, TT &&t_in, UT &&U_in) Construct a Hubbard Hamiltonian on a lattice.
Extends the Hückel model with on-site electron-electron repulsion: H = H_huckel + U sum_i n_{i,up} n_{i,down}
- Template Parameters:
EpsT – double or Eigen::VectorXd
TT – double or Eigen::MatrixXd
UT – double or Eigen::VectorXd
- Parameters:
lattice – Symmetric lattice graph defining the connectivity.
epsilon_in – On-site orbital energies. Scalar or VectorXd of size n.
t_in – Hopping integrals. Scalar or n x n MatrixXd.
U_in – On-site Coulomb repulsion. Scalar or VectorXd of size n.
- Throws:
std::invalid_argument – if U vector size mismatches the number of sites.
- Returns:
Tuple of (sparse one-body matrix, two-body map).
-
template<typename EpsT, typename TT>
inline Eigen::SparseMatrix<double> _build_huckel_integrals(const qdk::chemistry::data::LatticeGraph &lattice, EpsT &&epsilon_in, TT &&t_in) Construct a Hückel Hamiltonian on a lattice.
Builds the one-body Hamiltonian: H = sum_i epsilon_i n_i - sum_{<i,j>} t_ij (a_i^dag a_j + a_j^dag a_i) where n_i = sum_sigma a_{i,sigma}^dag a_{i,sigma} and the sum over <i,j> runs over edges of the lattice graph.
- Template Parameters:
EpsT – double or Eigen::VectorXd
TT – double or Eigen::MatrixXd
- Parameters:
lattice – Symmetric lattice graph defining the connectivity.
epsilon_in – On-site orbital energies. Scalar or VectorXd of size n.
t_in – Hopping integrals. Scalar or n x n MatrixXd.
- Throws:
std::invalid_argument – if dimensions mismatch, lattice is asymmetric, or empty.
- Returns:
Sparse one-body integral matrix (n x n).
-
template<typename EpsT, typename TT, typename UT, typename VT, typename ZT>
inline std::tuple<Eigen::SparseMatrix<double>, qdk::chemistry::data::SparseHamiltonianContainer::TwoBodyMap, double> _build_ppp_integrals(const qdk::chemistry::data::LatticeGraph &lattice, EpsT &&epsilon_in, TT &&t_in, UT &&U_in, VT &&V_in, ZT &&z_in) Construct a Pariser-Parr-Pople (PPP) Hamiltonian on a lattice.
Extends the Hubbard model with long-range intersite Coulomb interactions: H = H_hubbard + 1/2 sum_{i!=j} V_ij (n_i - z_i)(n_j - z_j)
Note
The 1/2 prefactor from the PPP formula is not included in the stored two-body integrals.
- Template Parameters:
EpsT – double or Eigen::VectorXd
TT – double or Eigen::MatrixXd
UT – double or Eigen::VectorXd
VT – double or Eigen::MatrixXd
ZT – double or Eigen::VectorXd
- Parameters:
lattice – Symmetric lattice graph defining the connectivity.
epsilon_in – On-site orbital energies. Scalar or VectorXd of size n.
t_in – Hopping integrals. Scalar or n x n MatrixXd.
U_in – On-site Coulomb repulsion. Scalar or VectorXd of size n.
V_in – Intersite Coulomb interaction matrix. Scalar or n x n MatrixXd.
z_in – Effective core charges. Scalar or VectorXd of size n.
- Throws:
std::invalid_argument – if V or z dimensions mismatch.
- Returns:
Tuple of (sparse one-body matrix, two-body map, energy offset).
-
inline const Eigen::MatrixXd &to_pair_param(const Eigen::MatrixXd &m, const qdk::chemistry::data::LatticeGraph &lattice, const std::string &name = "parameter")
Convert a per-pair parameter to MatrixXd with validation.
For MatrixXd input, validates that both dimensions match the number of lattice sites and returns a reference to the original matrix. For double input, broadcasts to a constant n x n MatrixXd.
- Parameters:
m – Per-pair parameter matrix.
lattice – Lattice graph whose site count defines the expected size.
name – Parameter name used in error messages.
- Throws:
std::invalid_argument – if the matrix dimensions do not match.
- Returns:
Reference to the validated matrix.
-
inline Eigen::MatrixXd to_pair_param(double val, const qdk::chemistry::data::LatticeGraph &lattice, const std::string& = "parameter")
Convert a scalar per-pair parameter to a constant n x n MatrixXd.
-
inline const Eigen::VectorXd &to_site_param(const Eigen::VectorXd &v, const qdk::chemistry::data::LatticeGraph &lattice, const std::string &name = "parameter")
Convert a per-site parameter to VectorXd with validation.
For VectorXd input, validates that the vector length matches the number of lattice sites and returns a reference to the original vector. For double input, broadcasts to a constant VectorXd.
- Parameters:
v – Per-site parameter vector.
lattice – Lattice graph whose site count defines the expected size.
name – Parameter name used in error messages.
- Throws:
std::invalid_argument – if the vector size does not match.
- Returns:
Reference to the validated vector.
-
inline Eigen::VectorXd to_site_param(double val, const qdk::chemistry::data::LatticeGraph &lattice, const std::string& = "parameter")
Convert a scalar per-site parameter to a constant VectorXd.
-
template<typename EpsT, typename TT, typename UT>
-
template<typename EpsT, typename TT, typename UT>