Projector Monte Carlo / FCIQMC
The purpose of Projector Monte Carlo is to stochastically sample the ground state, i.e. the eigenvector corresponding to the lowest eigenvalue of a quantum Hamiltonian, or more generally, a very large matrix. Rimu implements a flavor of Projector Monte Carlo called Full Configuration Interaction Quantum Monte Carlo (FCIQMC).
ProjectorMonteCarloProblem
To run a projector Monte Carlo simulation you set up a problem with ProjectorMonteCarloProblem
and solve it with solve
. Alternatively you can init
it with to obtain a PMCSimulation
struct, step!
through time steps, and solve!
it to completion.
Rimu.ProjectorMonteCarloProblem
— TypeProjectorMonteCarloProblem(hamiltonian::AbstractHamiltonian; kwargs...)
Defines a problem to be solved by projector quantum Monte Carlo (QMC) methods, such as the the FCIQMC
algorithm.
Common keyword arguments and defaults:
time_step = 0.01
: Initial time step size.last_step = 100
: Controls the number of steps.target_walkers = 1_000
: Target for the 1-norm of the coefficient vector.start_at = starting_address(hamiltonian)
: Define the initial state vector(s). An $r × s$ matrix of state vectors can be passed where $r$ is the number of replicas and $s$ the number of spectral states. See alsodefault_starting_vector
.style = IsDynamicSemistochastic()
: TheStochasticStyle
of the simulation.initiator = false
: Whether to use initiators. Can betrue
,false
, or a validInitiatorRule
.threading
: Default is to use multithreading and/or MPI if available. Set totrue
to forcePDVec
for the starting vector,false
for serial computation; may be overridden bystart_at
.reporting_strategy = ReportDFAndInfo()
: How and when to report results, seeReportingStrategy
.post_step_strategy = ()
: Extract observables (e.g.ProjectedEnergy
), seePostStepStrategy
.n_replicas = 1
: Number of synchronised independent simulations.replica_strategy = NoStats(n_replicas)
: Which results to report from replica simulations, seeReplicaStrategy
.n_spectral = 1
: Number of targeted spectral states. Setn_spectral > 1
to find excited states.spectral_strategy = GramSchmidt(n_spectral)
: TheSpectralStrategy
used for orthogonalizing spectral states.
Example
julia> hamiltonian = HubbardReal1D(BoseFS(1,2,3));
julia> problem = ProjectorMonteCarloProblem(hamiltonian; target_walkers = 500, last_step = 100);
julia> simulation = solve(problem);
julia> simulation.success[]
true
julia> size(DataFrame(simulation))
(100, 9)
Further keyword arguments:
starting_step = 1
: Starting step of the simulation.wall_time = Inf
: Maximum time allowed for the simulation.simulation_plan = SimulationPlan(; starting_step, last_step, wall_time)
: Defines the duration of the simulation. Takes precedence overlast_step
andwall_time
.ζ = 0.08
: Damping parameter for the shift update.ξ = ζ^2/4
: Forcing parameter for the shift update.shift_strategy = DoubleLogUpdate(; target_walkers, ζ, ξ)
: How to update theshift
, seeShiftStrategy
.time_step_strategy = ConstantTimeStep()
: Adjust time step or not, seeTimeStepStrategy
.algorithm = FCIQMC(; shift_strategy, time_step_strategy)
: The algorithm to use. Currenlty onlyFCIQMC
is implemented.shift
: Initial shift value or collection of shift values. Determined by default from the Hamiltonian and the starting vectors.initial_shift_parameters
: Initial shift parameters or collection of initial shift parameters. Overridesshift
if provided.max_length = 2 * target_walkers + 100
: Maximum length of the vectors.display_name = "PMCSimulation"
: Name displayed in progress bar (viaProgressLogging
).metadata
: User-supplied metadata to be added to the report. Must be an iterable of pairs or aNamedTuple
, e.g.metadata = ("key1" => "value1", "key2" => "value2")
. All metadata is converted to strings.random_seed = true
: Provide and store a seed for the random number generator. If set totrue
, a new random seed is generated fromRandomDevice()
. If set to number, this number is used as the seed. This seed is used bysolve
(andinit
) to re-seed the default random number generator (consistently on each MPI rank) such thatsolve
ing the sameProjectorMonteCarloProblem
twice will yield identical results. If set tofalse
, no seed is used and consecutive random numbers are used.minimum_size = 2*num_spectral_states(spectral_strategy)
: The minimum size of the basis used to construct starting vectors for simulations of spectral states, ifstart_at
is not provided.
CommonSolve.init
— Functioninit(p::ExactDiagonalizationProblem, [algorithm]; kwargs...)
Initialize a solver for an ExactDiagonalizationProblem
p
with an optional algorithm
. Returns a solver instance that can be solved with solve
.
For a description of the keyword arguments, see the documentation for ExactDiagonalizationProblem
.
init(problem::ProjectorMonteCarloProblem; copy_vectors=true)::PMCSimulation
Initialise a Rimu.PMCSimulation
.
See also ProjectorMonteCarloProblem
, solve!
, solve
, step!
, Rimu.PMCSimulation
.
CommonSolve.solve
— Functionsolve(::ProjectorMonteCarloProblem)::PMCSimulation
Initialize and solve a ProjectorMonteCarloProblem
until the last step is completed or the wall time limit is reached.
See also init
, solve!
, step!
, Rimu.PMCSimulation
, and solve(::ExactDiagonalizationProblem)
.
CommonSolve.solve!
— Functionsolve!(sm::PMCSimulation; kwargs...)::PMCSimulation
Solve a Rimu.PMCSimulation
until the last step is completed or the wall time limit is reached.
To continue a previously completed simulation, set a new last_step
or wall_time
using the keyword arguments. Optionally, changes can be made to the replica_strategy
, the post_step_strategy
, or the reporting_strategy
.
Optional keyword arguments:
last_step = nothing
: Set the last step to a new value and continue the simulation.wall_time = nothing
: Set the allowed wall time to a new value and continue the simulation.reset_time = false
: Reset theelapsed_time
counter and continue the simulation.empty_report = false
: Empty the report before continuing the simulation.replica_strategy = nothing
: Change the replica strategy. Requires the number of replicas to match the number of replicas in the simulationsm
. Impliesempty_report = true
.post_step_strategy = nothing
: Change the post-step strategy. Impliesempty_report = true
.reporting_strategy = nothing
: Change the reporting strategy. Impliesempty_report = true
.metadata = nothing
: Add metadata to the report.
See also ProjectorMonteCarloProblem
, init
, solve
, step!
, Rimu.PMCSimulation
.
CommonSolve.step!
— Functionstep!(sm::PMCSimulation)::PMCSimulation
Advance the simulation by one step.
Calling solve!
will advance the simulation until the last step or the wall time is exceeded. When completing the simulation without calling solve!
, the simulation report needs to be finalised by calling Rimu.finalize_report!
.
See also ProjectorMonteCarloProblem
, init
, solve!
, solve
, Rimu.PMCSimulation
.
After solve
or solve!
have been called the returned PMCSimulation
contains the results of the projector Monte Carlo calculation.
PMCSimulation
and report as a DataFrame
Rimu.PMCSimulation
— TypePMCSimulation
Holds the state and the results of a projector quantum Monte Carlo (PMC) simulation. Is returned by init(::ProjectorMonteCarloProblem)
and solved with solve!(::PMCSimulation)
.
Obtain the results of a simulation sm
as a DataFrame with DataFrame(sm)
.
Fields
problem::ProjectorMonteCarloProblem
: The problem that was solvedstate::Rimu.ReplicaState
: The current state of the simulationreport::Rimu.Report
: The report of the simulationmodified::Bool
: Whether the simulation has been modifiedaborted::Bool
: Whether the simulation has been abortedsuccess::Bool
: Whether the simulation has been completed successfullymessage::String
: A message about the simulation statuselapsed_time::Float64
: The time elapsed during the simulation
See also state_vectors
, ProjectorMonteCarloProblem
, init
, solve!
.
The DataFrame
returned from DataFrame(::PMCSimulation)
contains the time series data from the projector Monte Carlo simulation that is of primary interest for analysis. Depending on the reporting_strategy
and other options passed as keyword arguments to ProjectorMonteCarloProblem
it can have different numbers of rows and columns. The rows correspond to the reported time steps (Monte Carlo steps). There is at least one column with the name :step
. Further columns are usually present with additional data reported from the simulation.
For the default option algorithm = FCIQMC(; shift_strategy, time_step_strategy)
with a single replica (n_replicas = 1
) and single spectral state, the fields :shift
, :norm
, :len
will be present as well as others depending on the style
argument and the post_step_strategy
.
If multiple replicas or spectral states are requested, then the relevant field names in the DataFrame
will have a suffix identifying the respective replica simulation, e.g. the shift
s will be reported as shift_r1s1
, shift_r2s1
, ...
Many tools for analysing the time series data obtained from a ProjectorMonteCarloProblem
are contained in the Module StatsTools
.