Namespace qdk::chemistry::utils
-
namespace utils
Enums
-
enum class LogLevel
Log level enumeration for QDK Chemistry logging.
Values:
-
enumerator trace
Most verbose logging.
-
enumerator debug
Debug information.
-
enumerator info
General information.
-
enumerator warn
Warning messages.
-
enumerator error
Error messages.
-
enumerator critical
Critical errors.
-
enumerator off
Disable logging.
-
enumerator trace
Functions
Computes the default number of valence orbitals and electrons for a given structure and wavefunction.
This function analyzes the provided molecular structure and wavefunction to determine the number of electrons and orbitals in valence shell
- Parameters:
wavefunction – Shared pointer to the wavefunction containing electronic structure information (structure is extracted from wavefunction->orbitals->basis_set)
charge – The total charge of the molecular system, which should be equal to the charge set in scf calculation
- Returns:
A pair where:
first: Number of valence electrons
second: Number of valence orbitals
-
template<typename T, typename Hasher = std::hash<T>>
inline std::size_t hash_combine(std::size_t seed, const T &v) Combines a hash value with the hash of another value.
This function implements a hash combination algorithm that takes an existing hash seed and combines it with the hash of a new value. It uses a modified version of the boost::hash_combine algorithm with the golden ratio constant for better hash distribution.
Note
The algorithm uses the magic constant 0x9e3779b9 (derived from the golden ratio) and bit shifting for good hash distribution properties.
- Template Parameters:
T – The type of value to hash.
Hasher – The hash function type (defaults to std::hash<T>).
- Parameters:
seed – The existing hash value to combine with.
v – The value to hash and combine.
- Returns:
The combined hash value.
-
template<typename T, typename ...Args>
inline std::size_t hash_combine(std::size_t seed, const T &v, Args&&... args) Variadic overload that combines a hash seed with multiple values.
This function recursively combines the hash of multiple values into a single hash value, processing values left-to-right.
Example:
std::size_t h = hash_combine(0, x, y, z); // Combines hashes of x, y, z
- Template Parameters:
T – The type of the first value to hash.
Args – The types of the remaining values to hash.
- Parameters:
seed – The existing hash value to combine with.
v – The first value to hash and combine.
args – The remaining values to hash and combine.
- Returns:
The combined hash value.
-
void log_trace_entering(const std::source_location &location = std::source_location::current())
Logs a standardized trace message for entering a function.
Automatically uses std::source_location to determine the calling function and logs a message in the form: “[context] Entering method_name”
This function is typically called via the QDK_LOG_TRACE_ENTERING() macro for convenience.
Example:
void compute() { QDK_LOG_TRACE_ENTERING(); // Logs "[qdk:chemistry:...:compute] Entering" // ... }
- Parameters:
location – Automatically provided by std::source_location::current()
Rotate molecular orbitals using a rotation vector.
This function takes Orbitals and applies orbital rotations using a rotation vector, which can be taken from stability analysis eigenvectors.
The rotation is performed by:
Unpacking the rotation vector into an anti-Hermitian matrix
Computing the unitary rotation matrix via matrix exponential
Applying the rotation to the molecular orbital coefficients
See also
Note
restricted_external can break spin symmetry and solve external instabilities of RHF/RKS.
Note
This function assumes aufbau filling for occupation numbers.
Note
Orbital energies are invalidated by rotation and set to null.
- Parameters:
orbitals – The Orbitals to rotate
rotation_vector – The rotation vector (eigenvector from stability analysis, corresponding to the lowest eigenvalue). See StabilityResult::eigenvector_format for detailed size and indexing requirements.
num_alpha_occupied_orbitals – Number of alpha occupied orbitals
num_beta_occupied_orbitals – Number of beta occupied orbitals
restricted_external – If true and orbitals are restricted, creates unrestricted orbitals with rotated coefficients for alpha spin and unrotated coefficients for beta spin. Default is false.
- Throws:
std::runtime_error – if rotation vector size is invalid
- Returns:
A new Orbitals object with rotated molecular orbital coefficients
-
inline std::string to_snake_case(const char *input)
Convert a PascalCase/camelCase string to snake_case at runtime.
This function inserts an underscore before each uppercase letter (except at position 0) and converts all letters to lowercase.
Examples:
”Ansatz” -> “ansatz”
”ConfigurationSet” -> “configuration_set”
”StabilityResult” -> “stability_result”
- Parameters:
input – Input string in PascalCase or camelCase
- Returns:
std::string containing the snake_case version
-
class ContextLogger
- #include <qdk/chemistry/utils/logger.hpp>
A logger wrapper that automatically prepends source context to messages.
This class wraps the global logger and automatically includes file/method context in every log message. It provides the same interface as spdlog::logger (trace, debug, info, warn, error, critical) but prepends source location info.
Use via the QDK_LOGGER() macro which captures the call site’s source location.
Example:
QDK_LOGGER().info("Message with {} args", 2); // Output: [timestamp] [info] [file:method] Message with 2 args
Public Functions
-
inline explicit ContextLogger(const std::source_location &loc)
-
inline void critical(const std::string &msg)
-
inline void debug(const std::string &msg)
-
inline void error(const std::string &msg)
-
inline void info(const std::string &msg)
-
inline void trace(const std::string &msg)
-
inline void warn(const std::string &msg)
-
inline explicit ContextLogger(const std::source_location &loc)
-
class Logger
- #include <qdk/chemistry/utils/logger.hpp>
Centralized logging utility wrapper around spdlog for QDK Chemistry.
This class provides a consistent interface for logging throughout the QDK Chemistry library, wrapping the spdlog library with project-specific defaults and conventions. It uses a single global logger instance for efficiency, while still providing per-file/function context in log messages.
Public Static Functions
-
static std::shared_ptr<spdlog::logger> get()
Get the global logger instance.
Returns the single global logger instance for QDK Chemistry. The logger is lazily initialized on first call and reused thereafter.
New logger is created with default settings:
Colored console output
Inherit global log level
Consistent timestamped format
- Returns:
Shared pointer to the global logger instance
-
static LogLevel get_global_level()
Get the current global log level.
Returns the current global logging level. This uses mutex protection to ensure thread safety.
- Returns:
The current global log level
-
static std::string get_source_context(const std::source_location &location = std::source_location::current())
Get a formatted source context string for the given location.
Returns a string like “qdk:chemistry:utils:logger:method_name” that identifies the source file and method. This is automatically used by the ContextLogger to prefix log messages with per-file context.
- Parameters:
location – Source location (defaults to caller’s location)
- Returns:
Formatted context string
-
static void set_global_level(LogLevel level)
Set the global log level for all loggers.
Changes the logging level for the global logger instance. Messages below this level will be suppressed.
- Parameters:
level – The minimum log level to output
-
static std::shared_ptr<spdlog::logger> get()
-
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>
-
enum class LogLevel