API
Fermion Determinant Matrix
SmoQyElPhQMC.FermionDetMatrix — Type
FermionDetMatrix{T<:Number, E<:AbstractFloat}A abstract type to represent fermion determinant matrix
\[M = \left(\begin{array}{ccccc} I & & & & B_{0}\\ -B_{1} & I\\ & -B_{2} & \ddots\\ & & \ddots & \ddots\\ & & & -B_{L_{\tau}-1} & I \end{array}\right),\]
where $B_l$ are propagator matrices for imaginary-time slice $\tau = \Delta\tau \cdot l$ given an inverse temperature $\beta = \Delta\tau \cdot L_\tau$. A Fermion determinant matrix $M$ will be $N L_\tau \times N L_\tau$, where each propagator matrix $B_l$ is $N \times N$, where $N$ is the number of orbitals in the lattice.
SmoQyElPhQMC.SymFermionDetMatrix — Type
SymFermionDetMatrix{T<:Number, E<:AbstractFloat} <: FermionDetMatrix{T,E}A type to represent fermion determinant matrix
\[M = \left(\begin{array}{ccccc} I & & & & B_{0}\\ -B_{1} & I\\ & -B_{2} & \ddots\\ & & \ddots & \ddots\\ & & & -B_{L_{\tau}-1} & I \end{array}\right),\]
where
\[B_l = \left[ e^{-\Delta\tau K_l/2} \right]^\dagger e^{-\Delta\tau V_l} e^{-\Delta\tau K_l/2}\]
are Hermitian (symmetric if real) propagator matrices for imaginary-time slice $\tau = \Delta\tau \cdot l$ given an inverse temperature $\beta = \Delta\tau \cdot L_\tau$. A Fermion determinant matrix $M$ will be $N L_\tau \times N L_\tau$, where each propagator matrix $B_l$ is $N \times N$, where $N$ is the number of orbitals in the lattice. Here the matrix $e^{-\Delta\tau K_l/2}$ is approximated using the non-hermitian checkerboard approximation.
SmoQyElPhQMC.SymFermionDetMatrix — Method
SymFermionDetMatrix(
fermion_path_integral::FermionPathIntegral{T, E};
maxiter::Int = (fermion_path_integral.N * fermion_path_integral.Lτ),
tol::E = 1e-6
) where {T<:Number, E<:AbstractFloat}Initialize an instance of the SymFermionDetMatrix type.
SmoQyElPhQMC.AsymFermionDetMatrix — Type
AsymFermionDetMatrix{T<:Number, E<:AbstractFloat} <: FermionDetMatrix{T, E}A type to represent fermion determinant matrix
\[M = \left(\begin{array}{ccccc} I & & & & B_{0}\\ -B_{1} & I\\ & -B_{2} & \ddots\\ & & \ddots & \ddots\\ & & & -B_{L_{\tau}-1} & I \end{array}\right),\]
where
\[B_l = e^{-\Delta\tau V_l} e^{-\Delta\tau K_l}\]
are Hermitian (symmetric if real) propagator matrices for imaginary-time slice $\tau = \Delta\tau \cdot l$ given an inverse temperature $\beta = \Delta\tau \cdot L_\tau$. A Fermion determinant matrix $M$ will be $N L_\tau \times N L_\tau$, where each propagator matrix $B_l$ is $N \times N$, where $N$ is the number of orbitals in the lattice. Note that $e^{-\Delta\tau K_l}$ is represented using the non-hermitian checkerboard approximation.
SmoQyElPhQMC.AsymFermionDetMatrix — Method
AsymFermionDetMatrix(
fermion_path_integral::FermionPathIntegral{T, E};
maxiter::Int = (fermion_path_integral.N * fermion_path_integral.Lτ),
tol::E = 1e-6
) where {T<:Number, E<:AbstractFloat}Initialize an instance of the AsymFermionDetMatrix type.
Preconditioners
SmoQyElPhQMC.KPMPreconditioner — Type
abstract type KPMPreconditioner{T<:Number, E<:AbstractFloat} endAbstract type representing the KPM preconditioner appearing in the left-preconditioned linear system
\[P^{-1} \cdot \left[ M^{\dagger} M^{\phantom\dagger} \right] \cdot x = P^{-1} \cdot b\]
that is solved for using the Conjugate Gradient method, where $M$ is the fermion determinant matrix. Here, a KPM preconditioner represents
\[P^{-1} = \left[ \bar{M}^{\dagger} \bar{M}^{\phantom\dagger} \right]^{-1}\]
with a Chebyshev expansion in powers of
\[\bar{B} = \frac{1}{L_\tau} \sum_{l=0}^{L_\tau-1} B_l,\]
where
\[\bar{M} = \left(\begin{array}{ccccc} I & & & & \bar{B}\\ -\bar{B} & I\\ & -\bar{B} & \ddots\\ & & \ddots & \ddots\\ & & & -\bar{B} & I \end{array}\right).\]
SmoQyElPhQMC.KPMPreconditioner — Method
KPMPreconditioner(
fermion_det_matrix::FermionDetMatrix{T};
# Keyword Arguments
rng::AbstractRNG = Random.default_rng(),
rbuf::E = 0.10,
n::Int = 20,
a1::E = 1.0,
a2::E = 1.0
) where {T<:Number, E<:AbstractFloat}Initialize and return an instance of either the SymKPMPreconditioner or AsymKPMPreconditioner type.
Arguments
fermion_det_matrix::FermionDetMatrix{T}: Fermion determinant matrix.
Keyword Arguments
rng::AbstractRNG = Random.default_rng(): Random number generator.rbuf::E = 0.10: Relative buffer applied to eigenvalue bounds of $\bar{B}$ calculated by Lanczos.n::Int = 20: Number of Lanczos iterations used to approximate eigenvalue bounds.a1::E = 1.0: Controls maximum order of kpm expansion.a2::E = 1.0: Controls minimum order of kpm expansion.
SmoQyElPhQMC.SymKPMPreconditioner — Type
mutable struct SymKPMPreconditioner{T, E, Tfft, Tifft} <: KPMPreconditioner{T, E}Type representing the KPM preconditioner appearing in the left-preconditioned linear system
\[P^{-1} \cdot \left[ M^{\dagger} M^{\phantom\dagger} \right] \cdot x = P^{-1} \cdot b\]
that is solved for using the Conjugate Gradient method, where $M$ is the fermion determinant matrix defined using the symmetric propagator definition
\[B_l = e^{-\Delta\tau K_l/2} e^{-\Delta\tau V_l} e^{-\Delta\tau K_l/2}.\]
Here, a KPM preconditioner represents
\[P^{-1} = \left[ \bar{M}^{\dagger} \bar{M}^{\phantom\dagger} \right]^{-1}\]
with a Chebyshev expansion in powers of
\[\bar{B} = \frac{1}{L_\tau} \sum_{l=0}^{L_\tau-1} B_l,\]
where
\[\bar{M} = \left(\begin{array}{ccccc} I & & & & \bar{B}\\ -\bar{B} & I\\ & -\bar{B} & \ddots\\ & & \ddots & \ddots\\ & & & -\bar{B} & I \end{array}\right).\]
SmoQyElPhQMC.AsymKPMPreconditioner — Type
mutable struct AsymKPMPreconditioner{T, E, Tfft, Tifft} <: KPMPreconditioner{T, E}Type representing the KPM preconditioner appearing in the left-preconditioned linear system
\[P^{-1} \cdot \left[ M^{\dagger} M^{\phantom\dagger} \right] \cdot x = P^{-1} \cdot b\]
that is solved for using the Conjugate Gradient method, where $M$ is the fermion determinant matrix defined using the asymmetric propagator definition
\[B_l = e^{-\Delta\tau V_l} e^{-\Delta\tau K_l}.\]
Here, a KPM preconditioner represents
\[P^{-1} = \left[ \bar{M}^{\dagger} \bar{M}^{\phantom\dagger} \right]^{-1}\]
with a Chebyshev expansion in powers of
\[\bar{B} = \frac{1}{L_\tau} \sum_{l=0}^{L_\tau-1} B_l,\]
where
\[\bar{M} = \left(\begin{array}{ccccc} I & & & & \bar{B}\\ -\bar{B} & I\\ & -\bar{B} & \ddots\\ & & \ddots & \ddots\\ & & & -\bar{B} & I \end{array}\right).\]
Monte Carlo Update Methods
SmoQyElPhQMC.PFFCalculator — Type
PFFCalculator{T<:AbstractFloat}The PFFCalculator type, short for pseudo-fermion field calculator, is for facilitating the sampling the pseudo-fermion fields $\Phi$, evaluate the fermionic action $S_f$ and calculating it's partial derivatives $\partial S_f/\partial x_{\tau,i}$ with respect to each phonon field $x_{\tau,i}.$
SmoQyElPhQMC.PFFCalculator — Method
PFFCalculator(
# Arguments
electron_phonon_parameters::ElectronPhononParameters{T},
fermion_det_matrix::FermionDetMatrix{T};
) where {T<:Number, E<:AbstractFloat}Initialize an instance of the PFFCalculator type used for calculating the pseudo-fermion fields $\Phi$. The tol and maxiter keywords specify the tolerance and maximum number of iterations used when performing conjugate gradient solves to evaluate the fermionic action given the current fields $\Phi$.
SmoQyElPhQMC.EFAPFFHMCUpdater — Type
struct EFAPFFHMCUpdater{T<:AbstractFloat, PFFT, PIFFT}Type to define how to perform an Hybrid/Hamiltonian Monte Carlo (HMC) updates of the phonon fields using a fermionic action formulated by introducing a complex pseudofermion field (PFF) \Phi. Exact Fourier Acceleration (EFA) is also used to more efficiently sample the phonon fields.
SmoQyElPhQMC.EFAPFFHMCUpdater — Method
EFAPFFHMCUpdater(;
# Keyword Arguments
electron_phonon_parameters::ElectronPhononParameters{T},
Nt::Int,
Δt::E = π/(2*Nt),
η::E = 0.0,
δ::E = 0.05
) where {T<:Number, E<:AbstractFloat}Initialize an instance of EFAPFFHMCUpdater type, defining an EFA-PFF-HMC update.
Keyword Arguments
electron_phonon_parameters::ElectronPhononParameters{T}: Parameters defining the electron-phonon model.Nt::Int: Number of HMC time-steps.Δt::E: Time-step for HMC update.η::E = 0.0: Regularization parameter for EFA.δ::E = 0.05: Fractional noise to add to the time-stepΔt.
SmoQyDQMC.hmc_update! — Function
hmc_update!(
# ARGUMENTS
electron_phonon_parameters::ElectronPhononParameters{T,E},
hmc_updater::EFAPFFHMCUpdater{E};
# KEYWORD ARGUMENTS
fermion_path_integral::FermionPathIntegral{T,E},
fermion_det_matrix::FermionDetMatrix{T,E},
pff_calculator::PFFCalculator{E},
rng::AbstractRNG,
recenter!::Function = identity,
Nt::Int = hmc_updater.Nt,
Δt::E = hmc_updater.Δt,
δ::E = hmc_updater.δ,
tol_action::E = fermion_det_matrix.cg.tol,
tol_force::E = sqrt(fermion_det_matrix.cg.tol),
max_iter::Int = fermion_det_matrix.cg.maxiter,
preconditioner = I
) where {T, E}Perform an EFA-PFF-HMC update to the phonon fields. Acronym EFA-PFF-HMC stands for pseudofermion field (PPF) Hamiltonian/hybrid Monte Carlo (HMC) update with exact Fourier acceleration (EFA) used to reduce autocorrelation times.
Keyword Arguments with Default Values
recenter!::Function = identity: Function to recenter the phonon fields after the update.Nt::Int = hmc_updater.Nt: Number of HMC time-steps.Δt::E = hmc_updater.Δt: Time-step for HMC update.δ::E = hmc_updater.δ: Fractional noise to add to the time-stepΔt.tol_action::E = fermion_det_matrix.cg.tol: Tolerance used in CG solve to evaluate fermionic action.tol_force::E = sqrt(fermion_det_matrix.cg.tol): Tolerance used in CG solve to evaluate derivative of fermionic action.max_iter::Int = fermion_det_matrix.cg.maxiter: Maximum number of iterations for CG solve.preconditioner = I: Preconditioner used in CG solves.
SmoQyDQMC.reflection_update! — Function
reflection_update!(
# ARGUMENTS
electron_phonon_parameters::ElectronPhononParameters{T,E},
pff_calculator::PFFCalculator{E};
# KEYWORD ARGUMENTS
fermion_path_integral::FermionPathIntegral{T,E},
fermion_det_matrix::FermionDetMatrix{T,E},
rng::AbstractRNG,
preconditioner = I,
tol::E = fermion_det_matrix.cg.tol,
maxiter::Int = fermion_det_matrix.cg.maxiter,
phonon_types = nothing
) where {T<:Number, E<:AbstractFloat}Randomly sample a phonon mode in the lattice, and propose an update that reflects all the phonon fields associated with that phonon mode $x \rightarrow -x.$ The argument phonon_types specifies the phonon ID's that are included for randomly sampling a phonon mode in the lattice to perform a swap update on. If phonon_types = nothing, then all types of phonon modes are included. This function returns a tuple containing (accepted, iters), where accepted is a boolean indicating whether the update was accepted or rejected, and iters is the number of CG iterations performed to calculate the fermionic action.
SmoQyDQMC.swap_update! — Function
swap_update!(
# ARGUMENTS
electron_phonon_parameters::ElectronPhononParameters{T,E},
pff_calculator::PFFCalculator{E};
# KEYWORD ARGUMENTS
fermion_path_integral::FermionPathIntegral{T,E},
fermion_det_matrix::FermionDetMatrix{T,E},
rng::AbstractRNG,
preconditioner = I,
tol::E = fermion_det_matrix.cg.tol,
maxiter::Int = fermion_det_matrix.cg.maxiter,
phonon_type_pairs = nothing
) where {T<:Number, E<:AbstractFloat}Randomly sample a pairs of phonon modes and exchange the phonon fields associated with the pair of phonon modes. The argument phonon_type_pairs specifies pairs phonon IDs that are used to randomly samples a pairs of phonon modes. If phonon_type_pairs = nothing, then all possible pairs of phonon types/IDs are allowed. This function returns a tuple containing (accepted, iters), where accepted is a boolean indicating whether the update was accepted or rejected, and iters is the number of CG iterations performed to calculate the fermionic action.
SmoQyDQMC.radial_update! — Function
radial_update!(
# ARGUMENTS
electron_phonon_parameters::ElectronPhononParameters{T,E},
pff_calculator::PFFCalculator{E};
# KEYWORD ARGUMENTS
fermion_path_integral::FermionPathIntegral{T,E},
fermion_det_matrix::FermionDetMatrix{T,E},
rng::AbstractRNG,
preconditioner = I,
tol::E = fermion_det_matrix.cg.tol,
maxiter::Int = fermion_det_matrix.cg.maxiter,
phonon_id::Union{Nothing,Int} = nothing,
σ::E = 1.0
) where {T<:Number, E<:AbstractFloat}Perform a radial update to the phonon fields, as described by Algorithm 1 in the paper arXiv:2411.18218. Specifically, the proposed update to the phonon fields $x$ is a rescaling such that $x \rightarrow e^{\gamma} x$ where $\gamma \sim N(0, \sigma/\sqrt{d})$ and $d$ is the number of phonon fields being updated.
Measurement Methods
SmoQyElPhQMC.GreensEstimator — Type
GreensEstimator{
T<:AbstractFloat, D, Dp1, Dp3,
Gfft<:AbstractFFTs.Plan, Gifft<:AbstractFFTs.Plan,
Cfft<:AbstractFFTs.Plan, Cifft<:AbstractFFTs.Plan,
}This type is used to compute stochastic estimates of the Green's function and other correlation functions.
SmoQyElPhQMC.GreensEstimator — Method
GreensEstimator(
# Arguments
fermion_det_matrix::FermionDetMatrix{T,E},
model_geometry::ModelGeometry{D,E};
# Keyword Arguments
Nrv::Int = 10,
preconditioner = I,
rng::AbstractRNG = Random.default_rng(),
maxiter::Int = fermion_det_matrix.cgs.maxiter,
tol::E = fermion_det_matrix.cgs.tol
) where {D, T<:Number, E<:AbstractFloat}Initialize an instance of the type GreensEstimator.
Arguments
fermion_det_matrix::FermionDetMatrix{T,E}: Fermion determinant matrix.model_geometry::ModelGeometry{D,E}: Defines model geometry.
Keyword Arguments
Nrv::Int = 10: Number of random vectors used to approximate Green's function.preconditioner = I: Preconditioner used to solve linear system.rng = Random.default_rng(): Random number generator.maxiter::Int = fermion_det_matrix.cgs.maxiter: Maximum number of iterations for linear solver.tol::E = fermion_det_matrix.cgs.tol: Tolerance for linear solver.
SmoQyDQMC.make_measurements! — Function
make_measurements!(
measurement_container::NamedTuple,
fermion_det_matrix::FermionDetMatrix{T,E},
greens_estimator::GreensEstimator{E,D};
# Keyword Arguments Start Here
model_geometry::ModelGeometry{D,E},
fermion_path_integral::FermionPathIntegral{T,E},
tight_binding_parameters::TightBindingParameters{T,E},
electron_phonon_parameters::ElectronPhononParameters{T,E},
preconditioner = I,
rng::AbstractRNG = Random.default_rng(),
tol::E = fermion_det_matrix.cg.tol,
maxiter::Int = fermion_det_matrix.cg.maxiter
) where {T<:Number, E<:AbstractFloat, D}Make all measurements.
Chemical Potential Tuning
SmoQyDQMC.update_chemical_potential! — Function
update_chemical_potential!(
# ARGUMENTS
fermion_det_matrix::FermionDetMatrix{T,E},
greens_estimator::GreensEstimator{E,D};
# KEYWORD ARGUMENTS
chemical_potential_tuner::MuTunerLogger{E,T},
tight_binding_parameters::TightBindingParameters{T,E},
fermion_path_integral::FermionPathIntegral{T,E},
preconditioner = I,
rng::AbstractRNG = Random.default_rng(),
update_greens_estimator::Bool = true,
tol::E = fermion_det_matrix.cgs.tol,
maxiter::Int = fermion_det_matrix.cgs.maxiter
) where {D, T<:Number, E<:AbstractFloat}Update the chemical potential $\mu$ in the simulation to approach the target density/filling. If update_greens_estimator = true, then greens_estimator is initialized to reflect the current state of the fermion_det_matrix.