Reference

GreenFunctionMonteCarlo.AbstractConfigType
AbstractConfig{T,N}

An abstract type that represents a configuration (i.e. spins, bosons on a lattice). This type is a subtype of AbstractArray{T,N}.

Parameters

  • T: The type of the elements in the configuration.
  • N: The number of dimensions of the configuration.

Interface

  • Base.parent(x::AbstractConfig): return the parent array
  • Base.copy(x::AbstractConfig): create a full copy of the configuration
  • apply!(x::AbstractConfig, move::Any): apply a move to the configuration
  • fulfills_constraints(x::AbstractConfig, HilbertSpace::AbstractHilbertSpace): check if the configuration satisfies the constraints of the Hilbert space
source
GreenFunctionMonteCarlo.AbstractDiagObservableType

Abstract supertype for observables which are diagonal in the computational basis and may be measured for free in GFMC. An observable must be a function that takes a configuration and returns an array.

Interface:

  • obs(::O) returns the buffer array for the output of the observable.
  • O(out,Conf): Writes the observable for the given configuration to the preallocated array out. Returns out.
  • Base.copy(O::O): Returns a copy of the observable.

If no preallocated array is given, the observable defaults to using the buffer array obs(O).

source
GreenFunctionMonteCarlo.AbstractGuidingFunctionType
AbstractGuidingFunction

An abstract type for all guiding function implementations in Green Function Monte Carlo. Subtypes of this correspond to implementations of concrete guiding functions, particularly variational wavefunctions. For optimal performance, it is recommended to implement the function logpsidiff! for the specific guiding function type.

Interface

  • logψ(x::AbstractArray): return the logarithm of the guiding function evaluated at the configuration x.
  • logψ(x::AbstractArray, H::AbstractHilbertSpace): return the logarithm of the guiding function evaluated at the configuration x in the specified HilbertSpace.
  • log_psi_diff(x::AbstractArray, dx::AbstractArray, logψ::AbstractGuidingFunction, Buffer::AbstractGuidingFunctionBuffer, Hilbert::AbstractHilbertSpace): return the logarithm of the ratio of the guiding function evaluated at the configuration x and x+dx in the specified HilbertSpace. Returns -Inf if the move is not applicable.
  • get_params(logψ::AbstractGuidingFunction): return the parameters of the guiding function as a linear Array. It is recommended to use RecursiveArrayTools.jl for this purpose.

Interface (optional)

  • allocate_GWF_buffers(logψ::AbstractGuidingFunction, NBuffers::Integer): allocate NBuffers instances of a buffer for the guiding function. Defaults to an array of EmptyGWFBuffer instances.
  • compute_GWF_buffer!(Buffer::AbstractGuidingFunctionBuffer, logψ::AbstractGuidingFunction, x): compute the full buffer for the guiding function at the configuration x.
  • pre_move_affect!(Buffer::AbstractGuidingFunctionBuffer, x, moves, logψ::AbstractGuidingFunction): perform any necessary operations before computing ratios ψx´_ψx.
  • post_move_affect!(Buffer::AbstractGuidingFunctionBuffer, x, dx, logψ::AbstractGuidingFunction): perform any necessary operations after the move is applied.
  • HDF5.h5write(file::AbstractString, name::AbstractString, logψ::AbstractGuidingFunction): write the guiding function to an HDF5 file.
source
GreenFunctionMonteCarlo.AbstractGuidingFunctionBufferType
AbstractGuidingFunctionBuffer

An abstract type that serves as a base for defining guiding function buffer structures.

Interface

  • setBuffer!(BA::AbstractGuidingFunctionBuffer, BB::AbstractGuidingFunctionBuffer): set the buffer BA to the values of the buffer BB.
source
GreenFunctionMonteCarlo.AbstractHilbertSpaceType
AbstractHilbertSpace

An abstract type representing a Hilbert space. This type serves as a base for defining various specific Hilbert spaces used in the project. Generally, a hilbert space should be defined by the number of sites, the number of local degrees of freedom (i.e. spin), and the constraint that the configurations must satisfy.

Interface

  • constraint(HilbertSpace::AbstractHilbertSpace): Return the constraint that the configurations in the Hilbert space must satisfy.
  • Base.size(HilbertSpace::AbstractHilbertSpace): Return the size of a config in the Hilbert space.

Interface (optional)

  • fulfills_constraint(x,HilbertSpace::AbstractHilbertSpace): Check if the given configuration x satisfies the constraint of the specified HilbertSpace. Defaults to constraint(HilbertSpace)(x).
source
GreenFunctionMonteCarlo.AbstractMoveType
AbstractMove

Abstract type representing a move, i.e. an operation that changes a configuration x to a new configuration x'.

Abstract type representing a move operation in the context of Green Function Monte Carlo simulations.

This abstract type serves as a base for defining various move operations that can be performed during the simulation process. Specific move types should inherit from this abstract type and implement the required functionality.

Interface

  • apply!(x::AbstractConfig, move::AbstractMove): Apply the move to the configuration x.
  • isapplicable(x::AbstractConfig, move::AbstractMove, HilbertSpace::AbstractHilbertSpace): Check if the move is applicable to the configuration x within the specified HilbertSpace.
  • affected_sites(move::AbstractMove): Return the sites affected by the move.

It is further necessary to implement either of the three methods:

  • move_dx(move::AbstractMove) Return the values of the move applied to the configuration x. i.e. x′ - x.
  • move_dx_before(move::AbstractMove,x::AbstractConfig: Uses x to compute the move values before applying the move to the configuration. Defaults to move_dx(move)
  • move_dx_after(move::AbstractMove,x::AbstractConfig): Uses x to compute the move values after applying the move to the configuration. Defaults to move_dx(move)
source
GreenFunctionMonteCarlo.AbstractObservableType

Abstract supertype for observables. An observable must be a function that takes a spin configuration and returns an array. To define a subtype of O of AbstractObservable, one must define the following functions:

  • obs(::O) returns the buffer array for the output of the observable.
  • O(out,Conf): Writes the observable for the given configuration to the preallocated array out. Returns out.
  • Base.copy(O::O): Returns a copy of the observable.

If no preallocated array is given, the observable defaults to using the buffer array obs(O).

source
GreenFunctionMonteCarlo.AbstractObserverType

Abstract supertype for observables in GFMC which are recorded during the run. Observables should always contain the table recording reconfiguration processes, the energies, and the total weights of each Markov step.

Interface:

  • saveObservables_before!(Observables,i,Walkers,propagator): Saves the observables for the given iteration i and walker ensemble Walkers.
  • saveObservables_after!(Observables,i,Walkers,propagator): Saves the observables for the given iteration i and walker ensemble Walkers after reconfiguration.
source
GreenFunctionMonteCarlo.AbstractOperatorType
AbstractOperator

An abstract type representing a general operator. This serves as a base type for defining various specific operators in the context of the Green Function Monte Carlo project.

source
GreenFunctionMonteCarlo.AbstractPropagatorType
AbstractPropagator

An abstract type representing a propagator. Examples include discrete or continuous time methods.

Interface:

  • propagateWalkers!(WE::AbstractWalkerEnsemble, H::AbstractSignFreeOperator, logψ::AbstractGuidingFunction, Hilbert::AbstractHilbertSpace, propagator::AbstractPropagator, parallelization::AbstractParallelizationScheme, RNG::Random.AbstractRNG = Random.default_rng()): Propagate the walkers in the ensemble using the specified moves and parameters.
source
GreenFunctionMonteCarlo.AbstractReconfigurationSchemeType
AbstractReconfigurationScheme

An abstract type that serves as a base for defining different reconfiguration strategies in the context of the Green Function Monte Carlo framework.

Interface

  • reconfigurateWalkers!(Walkers::AbstractWalkerEnsemble,reconfiguration::AbstractReconfigurationScheme,rng::Random.AbstractRNG)
  • get_reconfigurationList(reconfiguration::AbstractReconfigurationScheme)
source
GreenFunctionMonteCarlo.AbstractSignFreeOperatorType
AbstractSignFreeOperator <: AbstractOperator

An abstract type representing a sign-free operator in the context of Green Function Monte Carlo simulations.

Interface:

  • get_diagonal(O::AbstractSignFreeOperator): return the diagonal operator associated with the sign-free operator O
  • get_offdiagonal_elements(O::AbstractSignFreeOperator): return the weights associated with the off-diagonal operator O
source
GreenFunctionMonteCarlo.AbstractWalkerEnsembleType
AbstractWalkerEnsemble

An abstract type that represents an ensemble of configurations (i.e. spins, bosons on a lattice for each walker) in the context of the Green Function Monte Carlo project.

Interface (required)

  • getConfig(X::AbstractWalkerEnsemble,α): get the configuration of the α-th walker in the ensemble.
  • getMoveWeights(X::AbstractWalkerEnsemble,α): get the weights of the moves for the α-th walker in the ensemble.
  • getBuffer(X::AbstractWalkerEnsemble,α): get the buffer associated with the α-th walker in the ensemble.
  • getWalkerWeights(X::AbstractWalkerEnsemble): get the weights of the Walkers in the ensemble.

Interface (optional)

  • Base.eachindex(X::AbstractWalkerEnsemble): iterate over the indices of the ensemble.
  • getNWalkers(X::AbstractWalkerEnsemble): get the number of walkers in the ensemble.
  • getConfigs(X::AbstractWalkerEnsemble): get all configurations in the ensemble.
source
GreenFunctionMonteCarlo.BasicAccumulatorType
struct BasicAccumulator{T_high<:AbstractFloat} <: AbstractObserver

A structure that represents a basic accumulator for observing energy and weights during simulations. Instead of storing observables at each step, the expectation values are computed on the fly, reducing the storage requirements significantly. Note that it is advisable to use a good guess of the average weight in the propagator to reduce numerical precision loss.

Type Parameters

  • T_high: The floating-point type used for observables, constrained to subtypes of AbstractFloat.
source
GreenFunctionMonteCarlo.BasicObserverType
struct BasicObserver{DT<:AbstractFloat} <: AbstractObserver

Saves the energies, weights and reconfigurations of the walkers during the Monte Carlo simulation.

Type Parameters

  • T: Data type for the energies and weights.
source
GreenFunctionMonteCarlo.BosonConfigType
struct BosonConfig{T, N, Arr<:AbstractArray{T, N}} <: AbstractConfig{T, N}

Represents a configuration of bosonic particles in a given formalism. This structure is parameterized by:

  • T: The type of the elements in the configuration (e.g., Bool or UInt8).
  • N: The dimensionality of the configuration space.
  • Arr<:AbstractArray{T, N}: The type of the array used to store the configuration data, which must be a subtype of AbstractArray with element type T and dimensionality N.

This type is a subtype of AbstractConfig{T, N}, which provides a common interface for configurations in the system and implements the AbstractArray interface.

source
GreenFunctionMonteCarlo.BosonHilbertSpaceType
struct BosonHilbertSpace{Cons<:AbstractConstraint} <: AbstractHilbertSpace

A structure representing the Hilbert space for bosonic systems.

Type Parameters

  • Cons<:AbstractConstraint: A type that defines the constraints applied to the bosonic Hilbert space. This allows for flexible customization of the space's properties.

Supertype

  • AbstractHilbertSpace: This structure is a subtype of AbstractHilbertSpace, indicating that it represents a specific type of quantum mechanical Hilbert space.
source
GreenFunctionMonteCarlo.CombinedObserverType
CombinedObserver{T}

A struct that represents a combined observer, which is a collection of multiple observers grouped together. This allows for observing multiple quantities or behaviors simultaneously.

Fields

  • Observers::T: A collection containing the individual observers.

Usage

This struct is useful for combining multiple observer objects into a single entity, enabling them to be managed and used collectively.

source
GreenFunctionMonteCarlo.ConfigObserverFunction
ConfigObserver(filename, config::AbstractConfig{T,N}, NSteps::Integer, NWalkers::Integer) where {T,N}

Creates a combined observer that tracks energy, weight, reconfigurations and configurations of the walkers during the Monte Carlo simulation. The observer saves the data to a file with the given filename.

Arguments

  • filename::String: The name of the file where the configuration data will be saved.
  • config::AbstractConfig{T,N}: The initial configuration object representing the system state.
  • NSteps::Integer: The number of simulation steps to be observed.
  • NWalkers::Integer: The number of walkers (or particles) in the simulation.
source
GreenFunctionMonteCarlo.ContinuousTimePropagatorType
ContinuousTimePropagator{T<:AbstractFloat} <: AbstractPropagator

A structure representing a imaginary time projection of a trial wavefunction, also referred to as the "Continuous time limit".

Type Parameters

  • T<:AbstractFloat: The floating-point type used for numerical computations (e.g., Float64, Float32).

Supertype

  • AbstractPropagator: This structure is a subtype of AbstractPropagator, indicating that it implements the required interface for propagators in the Green Function Monte Carlo framework.

Fields

  • dτ::T: The time step for the propagation, which is a floating-point value. A small value may be inefficient in exploring the Hilbert space, while a large value will lead to a more unstable propagation. A good starting point is dτ = 0.1.
  • w_avg_estimate::T: An estimate of the average weight to reduce floating point errors. Ideally given by the exact ground state energy of the system.
source
GreenFunctionMonteCarlo.DiagonalOperatorType
DiagonalOperator

An abstract type representing a diagonal operator in the context of Green Function Monte Carlo simulations. A diagonal operator is special in the sense that it will not change the configuration of the system when applied to it and will only return a number.

Interface:

  • (D::DiagonalOperator)(x) : return the value of the operator applied to the configuration x
  • (D::DiagonalOperator)(x, params): return the value of the operator applied to the configuration x with parameters params, which can be used to store buffers.
source
GreenFunctionMonteCarlo.EqualWeightSuperpositionType
struct EqualWeightSuperposition <: AbstractGuidingFunction

Represents a guiding function that models an equal weight superposition of states, i.e. ψ(x) = 1 for all x, provided that the configuration x is a valid configuration satisfying all constraints. This structure is a subtype of AbstractGuidingFunction and is used in the context of variational calculations within the Green Function Monte Carlo framework.

See Also

  • AbstractGuidingFunction: The abstract type that this struct extends.
source
GreenFunctionMonteCarlo.GFMCProblemType
struct GFMCProblem{WE<:AbstractWalkerEnsemble, Prop<:AbstractPropagator, GF<:AbstractGuidingFunction, 
                 SFO<:AbstractSignFreeOperator, HS<:AbstractHilbertSpace, PS<:AbstractParallelizationScheme, 
                 RS<:AbstractReconfigurationScheme} <: AbstractGFMCProblem

Represents a Green's Function Monte Carlo (GFMC) problem with the following type parameters:

  • WE<:AbstractWalkerEnsemble: The type of the walker ensemble used in the simulation.
  • Prop<:AbstractPropagator: The propagator responsible for evolving the system.
  • GF<:AbstractGuidingFunction: The guiding function used to improve the efficiency of the simulation.
  • SFO<:AbstractSignFreeOperator: The operator ensuring sign-free sampling in the simulation, containing only negative elements on the off-diagonal.
  • HS<:AbstractHilbertSpace: The Hilbert space in which the problem is defined.
  • PS<:AbstractParallelizationScheme: The parallelization scheme used for distributing computation.
  • RS<:AbstractReconfigurationScheme: The reconfiguration scheme used to manage walker populations.
source
GreenFunctionMonteCarlo.GreenFunctionMonteCarloModule

Overview

GreenFunctionMonteCarlo.jl is a Julia package designed for performing Green Function Monte Carlo (GFMC) simulations on lattice models.

Presently, this package treats only Hamiltonians that are free of the sign problem, i.e. whose elements satisfy

\[ H_{x, x'} \leq 0 \quad \forall x \neq x'\]

where $H_{x, x'}$ is the matrix element of the Hamiltonian between two configurations (spins or bosons) $x$ and $x'$.

Usage:

using GreenFunctionMonteCarlo, LinearAlgebra
NSites = 3
Nwalkers = 10
NSteps = 10
Hilbert = BosonHilbertSpace(NSites, HardCoreConstraint())
moves = eachcol(Bool.(I(NSites))) # each move flips a single spin
offdiagElements = -ones(NSites)
H = localOperator(eachrow(moves), offdiagElements, DiagOperator(x->0), Hilbert)

problem = GFMCProblem(BosonConfig(Hilbert), Nwalkers, ContinuousTimePropagator(0.1); logψ = EqualWeightSuperposition(), H, Hilbert)
Observer = ConfigObserver(BosonConfig(Hilbert), NSteps, Nwalkers) # Observer to measure the energy and configurations
runGFMC!(problem, NoObserver(), 100) #run for 100 steps without observing to equilibrate
runGFMC!(problem, Observer, NSteps) #run for NSteps steps
source
GreenFunctionMonteCarlo.HardCoreConstraintType
HardCoreConstraint <: AbstractConstraint

A struct representing a hard-core constraint in the system. This constraint enforces that certain configurations are not allowed, typically used in bosonic systems to model hard-core interactions where particles cannot occupy the same state or position.

This type is a subtype of AbstractConstraint, which serves as a base for defining various constraints in the system. Useful for simulating spin-1/2 systems.

source
GreenFunctionMonteCarlo.InverseMoveType
struct InverseMove{T<:AbstractMove} <: AbstractMove

A type representing an inverse move in a Monte Carlo simulation. This type is parameterized by T, which must be a subtype of AbstractMove. The InverseMove struct is used to encapsulate the concept of an inverse operation corresponding to a given move in the simulation.

Parameters

  • T: A subtype of AbstractMove that specifies the type of move for which this struct represents the inverse.
source
GreenFunctionMonteCarlo.JastrowType
struct Jastrow{T<:Real} <: AbstractGuidingFunction

A structure representing the Jastrow factor in a variational Monte Carlo simulation. This function strikes a good balance between accuracy and computational efficiency in the context of Green Function Monte Carlo. $log(ψ(x)) = \sum_i m_i x_i + \frac{1}{2}\sum_{i,j} v_{ij} x_i x_j$

Type Parameters

  • T<:Real: The numeric type used for the parameters of the Jastrow factor (e.g., Float64, Float32).

Usage

logψ = Jastrow(N,Float64) # Create a Jastrow factor for N sites with parameters of type Float64. logψ = Jastrow(conf,Float64) # Create a Jastrow factor for configurations similar to conf.

source
GreenFunctionMonteCarlo.NoObserverType
NoObserver <: AbstractObserver

An observer type that does not save anything. Can be used to evolve the system in the most efficient way, i.e. for equilibration.

source
GreenFunctionMonteCarlo.ObservableAccumulatorType
ObservableAccumulator{ObsType<:AbstractObservable, T_high<:AbstractFloat, T_low<:Real}

An accumulator for observables in Monte Carlo simulations. This struct is designed to collect and process measurements of a given observable type during the simulation, supporting both high-precision (T_high) and lower-precision (T_low) data types.

Type Parameters

  • ObsType<:AbstractObservable: The type of observable being accumulated.
  • T_high<:AbstractFloat: The floating-point type used for high-precision accumulation (e.g., Float64).
  • T_low<:Real: The type used for lower-precision to speed-up calculations (typically Float32).

Description

ObservableAccumulator is typically used as part of the observer pattern in Monte Carlo simulations, where it collects measurements of observables at each step and provides methods for statistical analysis, such as computing averages and variances.

See Also

source
GreenFunctionMonteCarlo.ProblemEnsembleType
ProblemEnsemble{P<:AbstractGFMCProblem} <: AbstractGFMCProblem

Run multiple GFMC problems in parallel. Useful to estimate errors.

Usage example:

P = ProblemEnsemble([GFMCProblem1, GFMCProblem2, ...])

problems = ProblemEnsemble([GFMCProblem(startConfig, NWalkers, ContinuousTimePropagator(dtau); logψ, H, Hilbert) for _ in 1:10])

Observers = [ConfigObserver("output_$i.h5",startConfig, NSteps, NWalkers) for i in 1:10] #note that each observer must have its own file

runGFMC!(problems, NoObserver(),100) #equilibrate
runGFMC!(problems, Observers)
source
GreenFunctionMonteCarlo.ProgressBarLoggerType
ProgressBarLogger(; kwargs...)

Creates a progress bar logger with customizable update intervals and additional options.

Arguments

  • kwargs...: Additional keyword arguments to customize the behavior of the progress bar logger. See GreenFunctionMonteCarlo.ProgressMeter.Progress for more details.

Returns

A progress bar logger instance that can be used to track and display progress in a task.

Example

source
GreenFunctionMonteCarlo.SimpleLoggerType
struct SimpleLogger <: AbstractLogger

A simple logger implementation that inherits from AbstractLogger.

Fields

  • n_report::Int: The number of steps in the log report.

This logger can be used to control and manage logging behavior in a straightforward manner.

source
GreenFunctionMonteCarlo.allocate_GWF_buffersFunction
allocate_GWF_buffers(logψ::AbstractGuidingFunction, NBuffers::Integer, x)

Allocates a specified number of guiding wave function (GWF) buffers.

Arguments

  • logψ::AbstractGuidingFunction: The guiding function for which the buffers are being allocated.
  • NBuffers::Integer: The number of buffers to allocate.
  • x: An exemplary configuration

Returns

  • An array filled with NBuffers instances of AbstractGuidingFunctionBuffer. Defaults to an array of NotImplementedBuffer instances.
source
GreenFunctionMonteCarlo.apply!Function
apply!(x::AbstractConfig, move::AbstractMove)

Throws a MethodError indicating that the apply! function has not been implemented for the given AbstractConfig type and move.

Arguments

  • x::AbstractConfig: An instance of a type that is a subtype of AbstractConfig.
  • move::AbstractMove: A move or operation to be applied to the configuration x.

Throws

  • MethodError: Always thrown to indicate that the method needs to be implemented for specific subtypes of AbstractConfig.
source
GreenFunctionMonteCarlo.continuos_time_propagation!Function
continuos_time_propagation!(WE::AbstractWalkerEnsemble, H::AbstractSignFreeOperator, logψ::AbstractGuidingFunction, Hilbert::AbstractHilbertSpace, dτ::Real, parallelization::MultiThreaded, RNG::Random.AbstractRNG = Random.default_rng())

Perform continuous time propagation on a walker ensemble for a fixed time step .

Arguments

  • WE::AbstractWalkerEnsemble: The ensemble of walkers to be propagated.
  • H::AbstractSignFreeOperator: The Hamiltonian operator used for propagation.
  • logψ::AbstractGuidingFunction: The guiding function for the propagation.
  • Hilbert::AbstractHilbertSpace: The Hilbert space in which the propagation occurs.
  • dτ::Real: The time step for the propagation.
  • w_avg_estimate::Real: An estimate of the average weight.
  • parallelization::MultiThreaded: Parallelization settings for the propagation.
  • RNG::Random.AbstractRNG: The random number generator to be used (default is Random.default_rng()).
source
GreenFunctionMonteCarlo.estimate_weights_continuousTime!Function
estimate_weights_continuousTime!(prob; Nepochs=5, Nsamples=100, mProj=50)

Estimate weights in a continuous-time Monte Carlo simulation. Useful for reducing the floating point errors for accumulators. The results may be passed to ContinuousTimePropagator(tau,w_avg_estimate) and BasicAccumulator(args...;weight_normalization)

Arguments

  • prob: The problem instance or data structure containing the simulation setup.
  • Nepochs: (Optional) Number of epochs to run the estimation. Default is 5.
  • Nsamples: (Optional) Number of samples per epoch. Default is 100.
  • mProj: (Optional) Number of projection steps or iterations. Default is 50.

Description

This function performs an in-place estimation of weights for a given problem using a continuous-time Monte Carlo approach. The process iterates over a specified number of epochs, drawing samples and projecting as configured by the keyword arguments.

Returns

  • Modifies prob in place which helps to equilibrate
  • Returns a rough estimate of the mean total weight and CT.wavgestimate.
source
GreenFunctionMonteCarlo.fulfills_constraintFunction
fulfills_constraint(x::AbstractArray, HilbertSpace::AbstractHilbertSpace)

Check if the given configuration x satisfies the constraint of the specified HilbertSpace.

Arguments

  • x::AbstractArray: The configuration to be checked.
  • HilbertSpace::AbstractHilbertSpace: The Hilbert space whose constraint need to be satisfied.

Returns

  • Bool: true if the configuration satisfies the constraint, false otherwise.
source
GreenFunctionMonteCarlo.getEnergiesFunction
getEnergies(weights, localEnergies, PMax; kwargs...)

Compute the energies based on the provided weights and local energies, for the projection steps p = 0,1,2,...,PMax.

Arguments

  • weights: A collection of weights associated with the samples.
  • localEnergies: A collection of local energy values corresponding to the samples.
  • PMax: An integer indicating the maximum projection order.

Keyword arguments

  • Gnp: optionally provided precomputed normalized accumulated weights.

Returns

Returns the computed energies based on the input parameters.

Notes

Ensure that the dimensions of weights and localEnergies are compatible.

source
GreenFunctionMonteCarlo.isapplicableFunction
isapplicable(x::AbstractConfig, move::AbstractArray, HilbertSpace::AbstractHilbertSpace)

Check if a given move is applicable to the current configuration within a specified Hilbert space.

Arguments

  • x::AbstractConfig: The current configuration.
  • move::AbstractMove: The proposed move to be applied.
  • HilbertSpace::AbstractHilbertSpace: The Hilbert space in which the configuration and move are defined.

Returns

  • Bool: true if the move is applicable, false otherwise.

Defaults to applying the move to the configuration, checking if the constraints are satisfied, and then reverting the move.

Make sure to implement this function for specific subtypes of AbstractConfig and AbstractHilbertSpace to ensure optimal performance. See also isapplicable(x::AbstractConfig, move::AbstractMove, constraint::AbstractConstraint).

source
isapplicable(x::AbstractConfig, move::AbstractMove, constraint::AbstractConstraint) -> Bool

Check if a given move is applicable to a configuration under a specified constraint.

Arguments

  • x::AbstractConfig: The configuration to which the move will be applied.
  • move::AbstractMove: The move that is being checked for applicability.
  • constraint::AbstractConstraint: The constraint that must be satisfied for the move to be applicable.

Returns

  • Bool: true if the move is applicable to the configuration under the given constraint, false otherwise.
source
GreenFunctionMonteCarlo.minimizeReconfiguration!Function
minimizeReconfiguration!(list)

given a list of reconfiguration indices, minimizes the number of reconfigurations by swapping elements in the list. Each walker that survives a reconfiguration step remains unchanged while walkers that are killed get assigned to a new index.

Arguments

  • list: A collection (e.g., an array) that will be reconfigured in-place.
source
GreenFunctionMonteCarlo.propagateWalkers!Function

Propagates walker ensemble X using a collection of moves according to the rules specified by the propagator P.

Usage:

propagateWalkers!(WE::AbstractWalkerEnsemble, H::AbstractSignFreeOperator, logψ::AbstractGuidingFunction, Hilbert::AbstractHilbertSpace, propagator::AbstractPropagator, parallelization::AbstractParallelizationScheme, RNG::Random.AbstractRNG = Random.default_rng())

source
GreenFunctionMonteCarlo.reconfigurateWalkers!Function

Performs an efficient reconfiguration of walkers. This reconfiguration will not remove walkers if they all have the same weight, which increases the efficiency as more walkers can contribute to the average.

Matteo Calandra Buonaura and Sandro Sorella Phys. Rev. B 57, 11446 (1998)

Arguments

  • Walkers::AbstractWalkerEnsemble: The ensemble of walkers to be reconfigured.
  • reconfigurationList: A list of indices that will be reconfigured.
  • rng::Random.AbstractRNG: The random number generator to be used.
source
GreenFunctionMonteCarlo.runGFMC!Function
runGFMC!(Walkers::AbstractWalkerEnsemble, Observables::AbstractObserver, reconfiguration::AbstractReconfigurationScheme, 
         range, propagator::AbstractPropagator, logψ::AbstractGuidingFunction, H::AbstractSignFreeOperator, 
         Hilbert::AbstractHilbertSpace, parallelizer::AbstractParallelizationScheme, logger::AbstractLogger, RNG::Random.AbstractRNG)

Runs the Green Function Monte Carlo (GFMC) simulation for a given system.

Arguments

  • Walkers::AbstractWalkerEnsemble: The ensemble of walkers representing the quantum state.
  • Observables::AbstractObserver: The observer object used to measure physical quantities during the simulation.
  • reconfiguration::AbstractReconfigurationScheme: The scheme used to reconfigure the walker ensemble to maintain population control.
  • range: The range of iterations or time steps for the simulation.
  • propagator::AbstractPropagator: The propagator used to evolve the walkers.
  • logψ::AbstractGuidingFunction: The guiding function (logarithmic form) used for importance sampling.
  • H::AbstractSignFreeOperator: The Hamiltonian operator of the system, assumed to be sign-free.
  • Hilbert::AbstractHilbertSpace: The Hilbert space on which the system is defined.
  • parallelizer::AbstractParallelizationScheme: The parallelization scheme used to distribute computation across resources.
  • logger::AbstractLogger: The logger object used to record simulation progress.
  • RNG::Random.AbstractRNG: The random number generator used for stochastic processes in the simulation.

Description

This function performs the GFMC simulation by evolving the walker ensemble using the provided propagator and guiding function. It measures observables diagonal in the computational Basis at each step and applies reconfiguration to control the walker population. The simulation is parallelized according to the specified parallelization scheme.

Notes

  • The function modifies the Walkers object in place.
  • Ensure that all input objects are properly initialized before calling this function.
source
runGFMC!(prob::GFMCProblem, Observables::AbstractObserver, range; 
         logger = default_logger(), rng = Random.default_rng(), 
         reconfiguration = prob.reconfiguration, Propagator = prob.Propagator, 
         logψ = prob.logψ, H = prob.H, Hilbert = prob.Hilbert, 
         parallelization = prob.parallelization)

Run the Green's Function Monte Carlo (GFMC) simulation for the given problem.

Arguments

  • prob::GFMCProblem: The GFMC problem instance containing the system configuration and parameters.
  • Observables::AbstractObserver: An observer object to track and record observables during the simulation.
  • range: Integer or range: The range of iterations or steps over which the simulation will be performed.

Keyword Arguments to override default values

  • logger: (Optional) A logger instance for logging simulation progress. Defaults to default_logger(). Use NoLogger() or nothing to disable logging.
  • rng: (Optional) A random number generator to ensure reproducibility. Defaults to Random.default_rng().
  • reconfiguration: (Optional) The reconfiguration method to be used during the simulation. Defaults to prob.reconfiguration.
  • Propagator: (Optional) The propagator function for the simulation. Defaults to prob.Propagator.
  • logψ: (Optional) The logarithm of the wavefunction. Defaults to prob.logψ.
  • H: (Optional) The Hamiltonian operator. Defaults to prob.H.
  • Hilbert: (Optional) The Hilbert space of the system. Defaults to prob.Hilbert.
  • parallelization: (Optional) The parallelization strategy for the simulation. Defaults to prob.parallelization.

Returns

This function modifies the prob and Observables in place to reflect the results of the simulation.

Notes

  • Ensure that the prob and Observables are properly initialized before calling this function.
  • If range is provided as an integer, it will be converted to a range 1:range.
  • If logger is nothing, it will default to NoLogger().
source