API

Rimu

Rimu.RimuModule
Rimu

Random Integrator for Many-Body Quantum Systems

source
Rimu.AllOverlapsType
AllOverlaps(n=2, operator=nothing) <: ReplicaStrategy{n}

Run n replicas and report overlaps between all pairs of replica vectors. If operator is not nothing, the overlap dot(c1, operator, c2) is reported as well. If operator is a tuple of operators, the overlaps are computed for all operators.

Column names in the report are of the form c{i}dotc{j} for vector-vector overlaps, and c{i}Op{k}c{j} for operator overlaps.

See ReplicaStrategy and AbstractHamiltonian (for an interface for implementing operators).

source
Rimu.DeltaMemoryType
DeltaMemory(Δ::Int) <: MemoryStrategy

Before updating the shift, memory noise with a memory length of Δ is applied, where Δ = 1 means no memory noise.

r̃ = (pnorm - tnorm)/(dτ*pnorm) + shift
r = r̃ - <r̃>
source
Rimu.DeltaMemory2Type
DeltaMemory2(Δ::Int) <: MemoryStrategy

Before updating the shift, memory noise with a memory length of Δ is applied, where Δ = 1 means no memory noise.

r̃ = pnorm - tnorm + shift*dτ*pnorm
r = (r̃ - <r̃>)/(dτ*pnorm)

The long-term average of r is not guaranteed to be zero.

source
Rimu.DeltaMemory3Type
DeltaMemory3(Δ::Int, level::Float64) <: MemoryStrategy

Before updating the shift, apply multiplicative memory noise with a memory length of Δ at level level, where Δ = 1 means no memory noise.

r̃ = (pnorm - tnorm)/pnorm + dτ*shift
r = r̃ - <r̃>
w .*= 1 + level*r
source
Rimu.DontUpdateType
DontUpdate(; targetwalkers = 1_000_000) <: ShiftStrategy

Don't update the shift. Return when targetwalkers is reached.

source
Rimu.DoubleLogProjectedType
DoubleLogProjected(; target, projector, ζ = 0.08, ξ = ζ^2/4) <: ShiftStrategy

Strategy for updating the shift according to the log formula with damping parameter ζ and ξ after projecting onto projector.

\[S^{n+1} = S^n -\frac{ζ}{dτ}\ln\left(\frac{P⋅Ψ^{(n+1)}}{P⋅Ψ^{(n)}}\right)-\frac{ξ}{dτ}\ln\left(\frac{P⋅Ψ^{(n+1)}}{\text{target}}\right)\]

source
Rimu.DoubleLogSumUpdateType
DoubleLogSumUpdate(; targetwalkers = 1000, ζ = 0.08, ξ = ζ^2/4, α = 1/2) <: ShiftStrategy

Strategy for updating the shift according to the log formula with damping parameters ζ and ξ.

\[S^{n+1} = S^n -\frac{ζ}{dτ}\ln\left(\frac{N_\mathrm{w}^{n+1}}{N_\mathrm{w}^n}\right) - \frac{ξ}{dτ}\ln\left(\frac{N_\mathrm{w}^{n+1}}{N_\mathrm{w}^\text{target}}\right),\]

where $N_\mathrm{w} =$ (1-α)*walkernumber() + α*UniformProjector()⋅ψ computed with walkernumber() and UniformProjector(). When ξ = ζ^2/4 this corresponds to critical damping with a damping time scale T = 2/ζ.

source
Rimu.DoubleLogUpdateType
DoubleLogUpdate(; targetwalkers = 1000, ζ = 0.08, ξ = ζ^2/4) <: ShiftStrategy

Strategy for updating the shift according to the log formula with damping parameter ζ and ξ.

\[S^{n+1} = S^n -\frac{ζ}{dτ}\ln\left(\frac{\|Ψ\|_1^{n+1}}{\|Ψ\|_1^n}\right)-\frac{ξ}{dτ}\ln\left(\frac{\|Ψ\|_1^{n+1}}{\|Ψ\|_1^\text{target}}\right)\]

When ξ = ζ^2/4 this corresponds to critical damping with a damping time scale T = 2/ζ.

source
Rimu.DoubleLogUpdateAfterTargetWalkersType
DoubleLogUpdateAfterTargetWalkers(targetwalkers, ζ = 0.08, ξ = 0.0016) <: ShiftStrategy

Strategy for updating the shift: After targetwalkers is reached, update the shift according to the log formula with damping parameter ζ and ξ. See DoubleLogUpdate.

source
Rimu.LogUpdateType
LogUpdate(ζ = 0.08) <: ShiftStrategy

Strategy for updating the shift according to the log formula with damping parameter ζ.

\[S^{n+1} = S^n -\frac{ζ}{dτ}\ln\left(\frac{\|Ψ\|_1^{n+1}}{\|Ψ\|_1^n}\right)\]

source
Rimu.LogUpdateAfterTargetWalkersType
LogUpdateAfterTargetWalkers(targetwalkers, ζ = 0.08) <: ShiftStrategy

Strategy for updating the shift: After targetwalkers is reached, update the shift according to the log formula with damping parameter ζ. See LogUpdate.

source
Rimu.MultiScalarType
MultiScalar

Wrapper over a tuple that supports +, -, min, and max. Used with MPI communication because SVectors are treated as arrays by MPI.Allreduce and Tuples do not support scalar operations.

Example

Suppose you want to compute the sum of a vector dv and also get the number of positive elements it has in a single pass. You can use MultiScalar:

julia> dv = DVec(:a => 1, :b => -2, :c => 1);

julia> s, p = mapreduce(+, values(dv)) do v
    Rimu.MultiScalar(v, Int(sign(v) == 1))
end;

julia> s, p
(0, 2)

This will work with MPIData.

Note that only MultiScalars with the same types can be operated on. This is a feature, as it forces type stability.

source
Rimu.ProjectedEnergyType
ProjectedEnergy(hamiltonian, projector; hproj=:vproj, vproj=:vproj) <: PostStepStrategy

After every step, compute hproj = dot(projector, hamiltonian, dv) and vproj = dot(projector, dv), where dv is the instantaneous coefficient vector. projector can be an AbstractDVec, or an AbstractProjector.

Reports to columns hproj and vproj, which can be used to compute projective energy, e.g. with Rimu.StatsTools.ratio_of_means. The keyword arguments hproj and vproj can be used to change the names of these columns. This can be used to make the names unique when computing projected energies with different projectors in the same run.

source
Rimu.ProjectedMemoryType
ProjectedMemory(Δ::Int, projector, pp::Number) <: MemoryStrategy
ProjectedMemory(Δ::Int, projector, v::AbstractDVec)

Before updating the shift, apply memory noise to minimize the fluctuations of the overlap of the coefficient vector with projector. Averaging over Δ time steps is applied, where Δ = 1 means no memory noise is applied. Use pp to initialise the value of the projection or pass v in order to initialise the projection with pp = projector.v.

r̃ = (projector⋅v - projector⋅w)/projector⋅v + dτ*shift
r = r̃ - <r̃>

where v is the coefficient vector before and w after applying a regular FCIQMC step.

source
Rimu.QMCStateType
QMCState

Holds all information needed to run FCIQMC, except the data frame. Holds a NTuple of ReplicaStates and various strategies that control the algorithm.

source
Rimu.ReplicaStateType
ReplicaState(v, w, pnorm, params, id)

Struct that holds all information needed for an independent run of the algorithm.

Can be advanced a step forward with advance!.

Fields

  • hamiltonian: the model Hamiltonian.
  • v: vector.
  • w: working memory.
  • pnorm: previous walker number (see walkernumber).
  • params: the FCIQMCRunStrategy.
  • id: appended to reported columns.
source
Rimu.ReplicaStrategyType
ReplicaStrategy{N}

An abstract type that controles how lomc! uses replicas. A subtype of ReplicaStrategy{N} operates on N replicas and must implement the following function:

Concrete implementations:

  • NoStats: run (possibly one) replica(s), but don't report any additional info.
  • AllOverlaps: report overlaps between all pairs of replica vectors.
source
Rimu.ReportDFAndInfoType
ReportDFAndInfo(; k=1, i=100, io=stdout, writeinfo=false) <: ReportingStrategy

The default ReportingStrategy. Report every kth step to a DataFrame and write info message to io every ith step (unless writeinfo == false). The flag writeinfo is useful for controlling info messages in MPI codes, e.g. by setting writeinfo=is_mpi_root().

source
Rimu.ReportToFileType
ReportToFile(; kwargs...) <: ReportingStrategy

Reporting strategy that writes the report directly to a file. Useful when dealing with long jobs or large numbers of replicas, when the report can incur a significant memory cost.

Keyword arguments

  • filename: the file to report to. If the file already exists, a new file is created.
  • chunk_size = 1000: the size of each chunk that is written to the file. A DataFrame of this size is collected in memory and written to disk. When saving, an info message is also printed to io.
  • save_if = true: if this value is true, save the report, otherwise ignore it. Use save_if=is_mpi_root() when running MPI jobs.
  • return_df: if this value is true, read the file and return the data frame at the end of computation. Otherwise, an empty DataFrame is returned.
  • io=stdout: The IO to print messages to. Set to devnull if you don't want to see messages printed out.
source
Rimu.RunTillLastStepType
RunTillLastStep(step::Int = 0 # number of current/starting timestep
             laststep::Int = 100 # number of final timestep
             shiftMode::Bool = false # whether to adjust shift
             shift = 0.0 # starting/current value of shift
             dτ::Float64 = 0.01 # current value of time step
) <: FciqmcRunStrategy

Parameters for running fciqmc!() for a fixed number of time steps. For alternative strategies, see FciqmcRunStrategy.

source
Rimu.ShiftMemoryType
ShiftMemory(Δ::Int) <: MemoryStrategy

Effectively replaces the fluctuating shift update procedure for the coefficient vector by an averaged shift over Δ timesteps, where Δ = 1 means no averaging.

source
Rimu.SignCoherenceType
SignCoherence(reference[; name=:coherence]) <: PostStepStrategy

After each step, compute the proportion of configurations that have the same sign as they do in the reference_dvec. Reports to a column named name, which defaults to coherence.

source
Rimu.TripleLogUpdateType
TripleLogUpdate(; targetwalkers = 1000, ζ = 0.08, ξ = ζ^2/4, η = 0.01) <: ShiftStrategy

Strategy for updating the shift according to the extended log formula with damping parameters ζ, ξ, and η.

\[S^{n+1} = S^n -\frac{ζ}{dτ}\ln\left(\frac{N_\mathrm{w}^{n+1}}{N_\mathrm{w}^n}\right) - \frac{ξ}{dτ}\ln\left(\frac{N_\mathrm{w}^{n+1}}{N_\mathrm{w}^\text{target}}\right) - \frac{η}{dτ}\ln\left(\frac{\|ℜ(Ψ^{n+1})\|_1^2 + \|ℑ(Ψ^{n+1})\|_1^2} {\|ℜ(Ψ^{n})\|_1^2 + \|ℑ(Ψ^{n})\|_1^2}\right),\]

where $N_\mathrm{w}$ is the walkernumber(). When ξ = ζ^2/4 this corresponds to critical damping with a damping time scale T = 2/ζ.

source
Rimu.WalkerLonelinessType
WalkerLoneliness(threshold=1) <: PostStepStrategy

After each step, compute the proportion of configurations that are occupied by at most threshold walkers. Reports to a column named loneliness.

source
Rimu._n_walkersMethod
_n_walkers(v, s_strat)

Returns an estimate of the expected number of walkers as an integer.

source
Rimu.advance!Method
advance!(report::Report, state::QMCState, replica::ReplicaState)

Advance the replica by one step. The state is used only to access the various strategies involved. Steps, stats, and computed quantities are written to the report.

Returns true if the step was successful and calculation should proceed, false when it should terminate.

source
Rimu.all_overlapsMethod
all_overlaps(operators, vectors)

Get all overlaps between vectors and operators. This function is overlpaded for MPIData.

source
Rimu.apply_memory_noise!Method
r = apply_memory_noise!(w, v, shift, dτ, pnorm, m_strat::MemoryStrategy)

Apply memory noise to w, i.e. w .+= r.*v, computing the noise r according to m_strat. Note that m_strat needs to be compatible with StochasticStyle(w). Otherwise, an error exception is thrown. See MemoryStrategy.

w is the walker array after fciqmc step, v the previous one, pnorm the norm of v, and r the instantaneously applied noise.

source
Rimu.fciqmc_col!Method
fciqmc_col!(w, ham, add, num, shift, dτ)
fciqmc_col!(::Type{T}, args...)
-> spawns, deaths, clones, antiparticles, annihilations

Spawning and diagonal step of FCIQMC for single column of ham. In essence it computes

w .+= (1 .+ dτ.*(shift .- ham[:,add])).*num.

Depending on T ==StochasticStyle(w), a stochastic or deterministic algorithm will be chosen. The possible values for T are:

source
Rimu.fciqmc_step!Function
fciqmc_step!(Ĥ, v, shift, dτ, pnorm, w;
                      m_strat::MemoryStrategy = NoMemory()) -> ṽ, w̃, stats

Perform a single matrix(/operator)-vector multiplication:

\[\tilde{v} = [1 - dτ(\hat{H} - S)]⋅v ,\]

where Ĥ == ham and S == shift. Whether the operation is performed in stochastic, semistochastic, or determistic way is controlled by the trait StochasticStyle(w). See StochasticStyle. w is a local data structure with the same size and type as v and used for working. Both v and w are modified.

Returns the result , a (possibly changed) reference to working memory , and the array stats = [spawns, deaths, clones, antiparticles, annihilations]. Stats will contain zeros when running in deterministic mode.

source
Rimu.finalize_report!Method
finalize_report!(::ReportingStrategy, report)

Finalize the report. This function is called after all steps in lomc! have finished.

source
Rimu.lomc!Method
lomc!(ham::AbstractHamiltonian, [v]; kwargs...) -> df, state
lomc!(state::QMCState, [df]; kwargs...) -> df, state

Linear operator Monte Carlo: Perform a projector quantum Monte Carlo simulation for determining the lowest eigenvalue of ham. v can be a single starting vector. The default choice is

v = DVec(starting_address(ham) => 10; style=IsStochasticInteger())

and triggers the integer walker FCIQMC algorithm. See DVec and StochasticStyle.

Keyword arguments, defaults, and precedence:

  • params::FciqmcRunStrategy = RunTillLastStep(laststep = 100, dτ = 0.01, shift = diagonal_element(ham, starting_address(ham))) - basic parameters of simulation state, see FciqmcRunStrategy; is mutated
  • laststep - can be used to override information otherwise contained in params
  • s_strat::ShiftStrategy = DoubleLogUpdate(targetwalkers = 100, ζ = 0.08, ξ = ζ^2/4) - how to update the shift, see ShiftStrategy
  • maxlength = 2 * s_strat.targetwalkers + 100 - upper limit on the length of v; when reached, lomc! will abort
  • style = IsStochasticInteger() - set StochasticStyle for default v; unused if v is specified.
  • post_step::NTuple{N,<:PostStepStrategy} = () - extract observables (e.g. ProjectedEnergy), see PostStepStrategy.
  • replica::ReplicaStrategy = NoStats(1) - run several synchronised simulation, see ReplicaStrategy.
  • r_strat::ReportingStrategy = ReportDFAndInfo() - how and when to report results, see ReportingStrategy
  • τ_strat::TimeStepStrategy = ConstantTimeStep() - adjust time step dynamically, see TimeStepStrategy
  • m_strat::MemoryStrategy = NoMemory() - experimental: inject memory noise, see MemoryStrategy
  • threading = :auto - can be used to control the use of multithreading (overridden by wm)
    • :auto - use multithreading if s_strat.targetwalkers ≥ 500
    • true - use multithreading if available (set shell variable JULIA_NUM_THREADS!)
    • false - run on single thread
  • wm - working memory; if set, it controls the use of multithreading and overrides threading; is mutated
  • df = DataFrame() - when called with AbstractHamiltonian argument, a DataFrame can be passed into lomc! that will be pushed into.

Return values

lomc! returns a named tuple with the following fields:

  • df: a DataFrame with all statistics being reported.
  • state: a QMCState that can be used for continuations.

Example

julia> add = BoseFS((1,2,3));


julia> hamiltonian = HubbardReal1D(add);


julia> df1, state = lomc!(hamiltonian);


julia> df2, _ = lomc!(state, df1; laststep=200); # Continuation run


julia> size(df1)
(100, 12)

julia> size(df2)
(200, 12)
source
Rimu.post_stepFunction
post_step(::PostStepStrategy, ::ReplicaState) -> kvpairs

Compute statistics after FCIQMC step. Should return a tuple of :key => value pairs. See also PostStepStrategy.

source
Rimu.refine_r_stratMethod
refine_r_strat(r_strat::ReportingStrategy) -> r_strat

Initialize the reporting strategy. This can be used to set up filenames or other attributes that need to be unique for a run of FCIQMC.

source
Rimu.replica_statsFunction
replica_stats(::ReplicaStrategy{N}, replicas::NTuple{N,ReplicaState}) -> (names, values)

Return the names and values of statistics reported by ReplicaStrategy. names should be a tuple of Symbols or Strings and values should be a tuple of the same length.

source
Rimu.report!Method
 report!(::ReportingStrategy, step, report::Report, keys, values, id="")
 report!(::ReportingStrategy, step, report::Report, nt, id="")

Report keys and values to report, which will be converted to a DataFrame before lomc! exits. Alternatively, a nt::NamedTuple can be passed in place of keys and values. If id is specified, it is appended to all keys. This is used to differentiate between values reported by different replicas.

To overload this function for a new ReportingStrategy, overload report!(::ReportingStrategy, step, args...) and apply the report by calling report!(args...).

source
Rimu.report!Method
report!(report, keys, values, id="")
report!(report, pairs, id="")

Write keys, values pairs to report that will be converted to a DataFrame later. Alternatively, a named tuple or a collection of pairs can be passed instead of keys and values.

The value of id is appended to the name of the column, e.g. report!(report, :key, value, :_1) will report value to a column named :key_1.

source
Rimu.report_after_stepMethod
report_after_step(::ReportingStrategy, step, report, state)

This function is called exactly once at the very end of a step. For example, it can be used to print some information to stdout.

source
Rimu.sort_into_targets!Method
sort_into_targets!(target, source, stats) -> agg, wm, agg_stats

Aggregate coefficients from source to agg and from stats to agg_stats according to thread- or MPI-level parallelism. wm passes back a reference to working memory.

source
Rimu.step_statsMethod
step_stats(::StochasticStyle)

Return a tuple of names (Symbol or String) and a zeros of values of the same length. These will be reported as columns in the DataFrame returned by lomc!.

source
Rimu.threshold_projected_deposit!Method
threshold_projected_deposit!

This function performs threshold projection before spawning, but only for IsDynamicSemistochastic with the project_later parameter set to false.

source
Rimu.update_dvec!Method
update_dvec!([::StochasticStyle,] dvec) -> dvec, nt

Perform an arbitrary transformation on dvec after the spawning step is completed and report statistics to the DataFrame.

Returns the new dvec and a NamedTuple nt of statistics to be reported.

When extending this function for a custom StochasticStyle, define a method for the two-argument call signature!

source
Rimu.update_dτMethod
update_dτ(s<:TimeStepStrategy, dτ, tnorm) -> new dτ

Update the time step according to the strategy s.

source
Rimu.update_shiftFunction
update_shift(s <: ShiftStrategy, shift, shiftMode, tnorm, pnorm, dτ, step, df, v_new, v_old)

Update the shift according to strategy s. See ShiftStrategy.

source

Reexported Submodules

Hamiltonians

Link to Module Hamiltionians.jl

BitStringAddresses

Link to Module BitStringAddresses.jl

DictVectors

Rimu.DictVectorsModule

Module that provides data structures that behave similar to sparse vectors, but are indexed by arbitrary types (could be non-integers) similarly to dictionaries. The idea is to do linear algebra with data structures that are neither subtyped to AbstractVector nor to AbstractDict and are suitable for use with KrylovKit.jl. For this, the abstract type and interface AbstractDVec is provided, with the concrete implementation of DVec

source
Rimu.DictVectors.AbstractDVecType
DictVectors.AbstractDVec{K,V}

Abstract type for sparse vectors with valtype V based on dictionary-like structures. The vectors are designed to work well with FCIQMC and KrylovKit.

They lie somewhere between AbstractDicts and sparse AbstractVectors, generally behaving like a dictionary, while supportting various linear algebra functionality. Indexing with a value not stored in the dictionary returns zero(V). Setting a stored value to 0 or below eps(V::AbstractFloat) removes the value from the dictionary. Their length signals the number of stored elements, not the size of the vector space.

They have a StochasticStyle which selects the spawning algorithm in FCIQMC.

To iterate over an AbstractDVec, use pairs or values.

Interface

The interface is similar to the AbstractDict interface.

Implement what would be needed for the AbstractDict interface (pairs, keys, values, setindex!, getindex, delete!, length, haskey, empty!, isempty) and, in addition:

source
Rimu.DictVectors.CoherentInitiatorType
CoherentInitiator(threshold) <: InitiatorRule

Initiator rule to be passed to InitiatorDVec. An initiator is a configuration add with a coefficient with magnitude abs(v[add]) > threshold. Rules:

  • Initiators can spawn anywhere.

  • Non-initiators can spawn to initiators.

  • Multiple non-initiators can spawn to a single non-initiator if their contributions add up to a value greater than the initiator threshold.

    See InitiatorRule.

source
Rimu.DictVectors.DVecType
DVec{K,V,D<:AbstractDict{K,V},S}

Dictionary-based vector-like data structure for use with FCIQMC and KrylovKit. While mostly behaving like a Dict, it supports various linear algebra operations such as norm and dot. It has a StochasticStyle that is used to select an appropriate spawning strategy in the FCIQMC algorithm.

See also: AbstractDVec.

Constructors

  • DVec(dict::AbstractDict[; style, capacity]): create a DVec with dict for storage. Note that the data may or may not be copied.

  • DVec(args...[; style, capacity]): args... are passed to the Dict constructor. The Dict is used for storage.

  • DVec{K,V}([; style, capacity]): create an empty DVec{K,V}.

  • DVec(dv::AbstractDVec[; style, capacity]): create a DVec with the same contents as adv. The style is inherited from dv by default.

The default style is selected based on the DVec's valtype (see default_style). If a style is given and the valtype does not match the style's eltype, the values are converted to an appropriate type.

The capacity argument is optional and sets the initial size of the DVec via sizehint!.

Examples

julia> dv = DVec(:a => 1)
DVec{Symbol,Int64} with 1 entries, style = IsStochasticInteger{Int64}()
  :a => 1

julia> dv = DVec(:a => 2, :b => 3; style=IsDynamicSemistochastic())
DVec{Symbol,Float64} with 2 entries, style = IsDynamicSemistochastic{Float64, true}(1.0, Inf, 1.0)
  :a => 2.0
  :b => 3.0
source
Rimu.DictVectors.InitiatorType
Initiator(threshold) <: InitiatorRule

Initiator rule to be passed to InitiatorDVec. An initiator is a configuration add with a coefficient with magnitude abs(v[add]) > threshold. Rules:

  • Initiators can spawn anywhere.
  • Non-initiators can spawn to initiators.

See InitiatorRule.

source
Rimu.DictVectors.InitiatorDVecType
InitiatorDVec{K,V} <: AbstractDVec{K,V}

Dictionary-based vector-like data structure for use with lomc! and KrylovKit.jl. See AbstractDVec. Functionally identical to DVec, but contains InitiatorValues internally in order to facilitate initiator methods. How the initiators are handled is controlled by the initiator keyword argument (see below).

Constructors

  • InitiatorDVec(dict::AbstractDict[; style, initiator, capacity]): create an InitiatorDVec with dict for storage. Note that the data may or may not be copied.

  • InitiatorDVec(args...[; style, initiator, capacity]): args... are passed to the Dict constructor. The Dict is used for storage.

  • InitiatorDVec{K,V}([; style, initiator, capacity]): create an empty InitiatorDVec{K,V}.

  • InitiatorDVec(dv::AbstractDVec[; style, initiator, capacity]): create an InitiatorDVec with the same contents as dv. The style is inherited from dv by default.

Keyword arguments

  • style: A valid StochasticStyle. The default is selected based on the InitiatorDVec's valtype (see default_style). If a style is given and the valtype does not match the style's eltype, the values are converted to an appropriate type.

  • initiator = Initiator(1): A valid InitiatorRule. See Initiator.

  • capacity: Indicative size as Int. Optional. Sets the initial size of the InitiatorDVec via sizehint!.

source
Rimu.DictVectors.IsDynamicSemistochasticType
IsDynamicSemistochastic{T=Float64}(rel_threshold=1, abs_threshold=Inf, proj_threshold=1) <: StochasticStyle{T}

QMC propagation with non-integer walker numbers and reduced noise. All possible spawns are performed deterministically when number of walkers in a configuration is high. Stochastic vector compression with threshold proj_threshold is applied after spawning and diagonal death steps.

Unlike with IsStochasticWithThreshold, when late_projection is set to true, walker annihilation is done before the stochastic vector compression.

Parameters:

  • late_projection = true: If set to true, threshold projection is done after all spawns are collected, otherwise, values are projected as they are being spawned.

  • rel_threshold = 1.0: If the walker number on a configuration times this threshold is greater than the number of offdiagonals, spawning is done deterministically. Should be set to 1 or more for best performance.

  • abs_threshold = Inf: If the walker number on a configuration is greater than this value, spawning is done deterministically. Can be set to e.g abs_threshold = 0.1 * target_walkers.

  • proj_threshold = 1.0: Values below this number are stochastically projected to this value or zero. See also IsStochasticWithThreshold.

See also StochasticStyle.

source
Rimu.DictVectors.IsStochastic2PopType
IsStochastic2Pop{T=Complex{Int}}() <: StochasticStyle{T}

Trait for generalised vector of configurations indicating stochastic propagation with complex walker numbers representing two populations of integer walkers.

See also StochasticStyle.

source
Rimu.DictVectors.IsStochasticWithThresholdType
IsStochasticWithThreshold(threshold=1.0) <: StochasticStyle

Trait for generalised vector of configurations indicating stochastic propagation with real walker numbers and cutoff threshold.

During stochastic propagation, walker numbers small than threshold will be stochastically projected to either zero or threshold.

See also StochasticStyle.

source
Rimu.DictVectors.SimpleInitiatorType
SimpleInitiator(threshold) <: InitiatorRule

Simplified initiator rule to be passed to InitiatorDVec. An initiator is a configuration add with a coefficient with magnitude abs(v[add]) > threshold. Rules:

  • Initiators can spawn anywhere.
  • Non-initiators cannot spawn.

See InitiatorRule.

source
Rimu.DictVectors.StochasticStyleType
StochasticStyle(v)

Abstract type. When called as a function it returns the native style of the generalised vector v that determines how simulations are to proceed.

Implemented styles

Usage

Concrete StochasticStyles can be used for the style keyword argument of lomc! and DVec.

Interface

When defining a new StochasticStyle, subtype it as MyStyle<:StochasticStyle{T} where T is the concrete value type the style is designed to work with.

For it to work with lomc!, a StochasticStyle must define the following:

Optionally, it can also define update_dvec!, which can be used to perform arbitrary transformations on the generalised vector after the spawning step is complete.

source
OrderedCollections.freezeMethod
freeze(dv)

Create a "frozen" version of dv which can no longer be modified or used in the conventional manner, but supports faster dot products.

If dv is an MPIData, synchronize its contents among the ranks first.

source
Rimu.DictVectors.deposit!Method
deposit!(w::AbstractDVec, add, val, parent::Pair)

Add val into w at address add, taking into account initiator rules if applicable. parent contains the address => value pair from which the pair add => val was created. InitiatorDVec can intercept this and add its own functionality.

source
Rimu.DictVectors.storageFunction
storage(dvec) -> AbstractDict

Return the raw storage associated with dvec as an AbstractDict. Used in MPI communication.

source
Rimu.DictVectors.walkernumberMethod
walkernumber(w)

Compute the number of walkers in w. It is used for updating the shift. Overload this function for modifying population control.

In most cases walkernumber(w) is identical to norm(w,1). For AbstractDVecs with complex coefficients it reports the one norm separately for the real and the imaginary part as a ComplexF64. See Norm1ProjectorPPop.

source
Rimu.DictVectors.zero!Method
zero!(v)

Replace v by a zero vector as an inplace operation. For AbstractDVec types it means removing all non-zero elements. For AbstractArrays, it sets all of the values to zero.

source

ConsistentRNG

Link to Module ConsistentRNG.jl

Unexported Submodules

StatsTools

Link to Module Rimu/StatsTools

Blocking

Rimu.BlockingModule

Blocking

Module that contains functions performing the Flyvbjerg-Petersen (J. Chem. Phys. 91, 461 (1989)) blocking analysis for evaluating the standard error on a correlated data set. A "M-test" is also implemented based on Jonsson (Phys. Rev. E 98, 043304, (2018)).

source
Rimu.Blocking.autoblockMethod
autoblock(df::DataFrame; start = 1, stop = size(df)[1])
-> s̄, σs, ē, σe, k

Determine mean shift and projected energy with respective standard errors σs and σe by blocking analsis from the DataFrame df returned from fciqmc!(). The number k of blocking steps and decorrelation time 2^k are obtained from the M-test for the shift and also applied to the projected energy, assuming that the projected quantities decorrelate on the same time scale. Only the real part of the shift is considered. Returns a named tuple.

source
Rimu.Blocking.autoblockMethod
autoblock(dftup::Tuple; start = 1, stop = size(dftup[1])[1])
-> s̄1, σs1, s̄2, σs2, ē1, σe1, ē2, σe2, ēH, σeH, k

Replica version. dftup is the tuple of DataFrames returned from replica fciqmc!(). Returns a named tuple with shifts and three variational energy estimators and respective errors obtained from blocking analysis. The larger of the k values from M-tests on the two shift time series is used.

source
Rimu.Blocking.autocovarianceMethod
autocovariance(v::Vector,h::Int; corrected::Bool=true)

$\hat{\gamma}(h) =\frac{1}{n}\sum_{t=1}^{n-h}(v_{t+h}-\bar{v})(v_t-\bar{v})^*$ Calculate the autocovariance of dataset v with a delay h. If corrected is true (the default) then the sum is scaled with n-h, whereas the sum is scaled with n if corrected is false where n = length(v).

source
Rimu.Blocking.blockAndMTestMethod
v̄, σ, σσ, k, df = blockAndMTest(v::Vector)

Perform a blocking analysis and M-test on v returning the mean , standard error σ, its error σσ, the number of blocking steps k, and the DataFrame df with blocking data.

source
Rimu.Blocking.blockerMethod
blocker(v::Vector) -> new_v::Vector

Reblock the data by successively taking the mean of two adjacent data points to form a new vector with a half of the length(v). The last data point will be discarded if length(v) is odd.

source
Rimu.Blocking.blockingMethod
blocking(x::Vector,y::Vector) -> df::DataFrame

Perform a blocking analysis for the quotient of means x̄/ȳ from two data sets. If corrected is true (the default) then the sums in both variance and covariance are scaled with n-1, whereas the sums are scaled with n if corrected is false where n = length(x) = length(y). Entries in returned dataframe:

  • blocks = number of blocks in current blocking step;
  • mean_x, SD_x, SE_x, SE_SE_x = the mean, standard deviation, standard error and error on standard error estimated for dataset x;
  • mean_y, SD_y, SE_y, SE_SE_y = ditto. for dataset y;
  • Covariance = the covariance between data in x and y;
  • mean_f = x̄/ȳ;
  • SE_f = standard error estimated for x̄/ȳ.
source
Rimu.Blocking.blockingMethod
blocking(v::Vector; corrected::Bool=true) -> df

Perform a blocking analysis according to Flyvberg and Peterson JCP (1989) for single data set and return a DataFrame with statistical data for each blocking step. M-test data according to Jonsson PRE (2018) is also provided. If corrected is true (the default) then the sum in var is scaled with n-1 and in autocovariance is scaled with n-h, whereas the sum is scaled with n for both if corrected is false where n = length(v).

source
Rimu.Blocking.blocking_oldMethod
blocking(v::Vector; typos = nothing) -> df

Perform a blocking analysis according to Flyvberg and Peterson JCP (1989) for single data set and return a DataFrame with statistical data for each blocking step. M-test data according to Jonsson PRE (2018) is also provided.

Keyword argument typos

  • typos = nothing - correct all presumed typos.
  • typos = :FP - use Flyvberg and Peterson (correct) standard error and Jonsson formul for M.
  • typos = :Jonsson - calculate M and standard error as written in Jonsson.
source
Rimu.Blocking.combination_divisionMethod
combination_division(x::Vector,y::Vector; corrected::Bool=true)

Find the standard error on the quotient of means x̄/ȳ from two data sets, note that the standard errors are different on $(x̄/ȳ) \neq \bar{(\frac{x}{y})}$. If corrected is true (the default) then the sums in both variance and covariance are scaled with n-1, whereas the sums are scaled with n if corrected is false where n = length(x) = length(y).

source
Rimu.Blocking.covarianceMethod
covariance(x::Vector,y::Vector; corrected::Bool=true)

Calculate the covariance between the two data sets x and y with equal length. If corrected is true (the default) then the sum is scaled with n-1, whereas the sum is scaled with n if corrected is false where n = length(x) = length(y).

source
Rimu.Blocking.gWMethod
gW(norm::AbstractArray, shift::AbstractArray, dt, [b]; pad = :true) -> g
gW(df::DataFrame, [b]; pad = :true) -> g

Compute the growth witness

\[G^{(n)} = S^{(n)} - \frac{\vert\mathbf{c}^{(n+1)}\vert - \vert\mathbf{c}^{(n)}\vert}{\vert\mathbf{c}^{(n)}\vert d\tau},\]

where S is the shift and $\vert\mathbf{c}^{(n)}\vert ==$ norm[n, 1]. Setting b ≥ 1 a sliding average over b time steps is computed.

If pad is set to :false then the returned array g has the length length(norm) - b. If set to :true then g will be padded up to the same length as norm and shift.

source
Rimu.Blocking.growthWitnessMethod
growthWitness(norm::AbstractArray, shift::AbstractArray, dt; b = 30, pad = :true) -> g
growthWitness(df::DataFrame; b = 30, pad = :true) -> g

Compute the growth witness

\[G_b^{(n)} = S̄^{(n)} - \frac{\log\vert\mathbf{c}^{(n+b)}\vert - \log\vert\mathbf{c}^{(n)}\vert}{b d\tau},\]

where is an average of the shift over b time steps and $\vert\mathbf{c}^{(n)}\vert ==$ norm[n]. The parameter b ≥ 1 averages the derivative quantity over b time steps and helps suppress noise.

If pad is set to :false then the returned array g has the length length(norm) - b. If set to :true then g will be padded up to the same length as norm and shift.

source
Rimu.Blocking.mtestMethod
mtest(df::DataFrame; warn = true) -> k

The "M test" based on Jonsson, M. Physical Review E, 98(4), 043304, (2018). Expects df to be output of a blocking analysis with column df.M containing relevant M_j values, which are compared to a χ^2 distribution. Returns the row number k where the M-test is passed. If the M-test has failed mtest() returns the value -1 and optionally prints a warning message.

source
Rimu.Blocking.seMethod
se(v::Vector;corrected::Bool=true)

Calculate the standard error of the dataset v. If corrected is true (the default) then the sum in std is scaled with n-1, whereas the sum is scaled with n if corrected is false where n = length(v).

source
Rimu.Blocking.smoothenMethod
smoothen(noisy::AbstractVector, b; pad = :true)

Smoothen the array noisy by averaging over a sliding window of length b. Pad to length(noisy) if pad == true. Otherwise, the returned array will have the length length(noisy) - b.

source

RMPI

Rimu.RMPI.MPIAllToAllType
 MPIAllToAll

All-to-all communication strategy. The communication works in two steps: first MPI.Alltoall! is used to communicate the number of walkers each rank wants to send to other ranks, then MPI.Alltoallv! is used to send the walkers around.

Constructor

  • MPIAllToAll(Type{P}, np, id, comm): Construct an instance with pair type P on np processes with current rank id.
source
Rimu.RMPI.MPIDataType
MPIData(data; kwargs...)

Wrapper used for signaling that this data is part of a distributed data structure and communication should happen with MPI.

Keyword arguments:

source
Rimu.RMPI.MPIDataIteratorType
MPIDataIterator{I,M<:MPIData}

Iterator over keys, values, or pairs of a dv::MPIData. Unlike its name would suggest, it does not actually support iteration. To perform computations with it, use mapreduce, or its derivatives (sum, prod, reduce...), which will perform the reduction accross MPI ranks.

source
Rimu.RMPI.MPINoWalkerExchangeType
MPINoWalkerExchange(nprocs, my_rank, comm)

Strategy for for not exchanging walkers between ranks. Consequently there will be no cross-rank annihilations.

source
Rimu.RMPI.MPIOneSidedType
MPIOneSided(nprocs, myrank, comm, ::Type{T}, capacity)

Communication buffer for use with MPI one-sided communication (remote memory access). Up to capacity elements of type T can be exchanged between MPI ranks via put. It is important that isbitstype(T) == true. Objects of type MPIOneSided have to be freed manually with a (blocking) call to free().

source
Rimu.RMPI.MPIPointToPointType
MPIPointToPoint{N,A}

Point-to-point communication strategy. Uses circular communication using MPI.Send and MPI.Recv!.

Constructor

  • MPIPointToPoint(::Type{P}, np, id, comm): Construct an instance with pair type P on np processes with current rank id.
source
Base.:*Method
*(lop::AbstractHamiltonian, md::MPIData)

Allocating "Matrix"-"vector" multiplication with MPI-distributed "vector" md. The result is similar to localpart(md) with all content having been communicated to the correct targetrank. MPI communicating.

See MPIData.

source
Base.lengthMethod
length(md::MPIData)

Compute the length of the distributed data on every MPI rank with MPI.Allreduce. MPI syncronizing.

source
LinearAlgebra.normFunction
norm(md::MPIData, p=2)

Compute the norm of the distributed data on every MPI rank with MPI.Allreduce. MPI syncronizing.

source
Rimu.ConsistentRNG.check_crng_independenceMethod
ConsistentRNGs.check_crng_independence(dv::MPIData)

Does a sanity check to detect dependence of random number generators across all MPI ranks. Returns the size of the combined RNG state, i.e. mpi_size()*Threads.nthreads()*fieldcount(ConsistentRNG.CRNG). MPI syncronizing.

source
Rimu.ConsistentRNG.sync_cRandnMethod
sync_cRandn(md::MPIData)

Generate one random number with cRandn() in a synchronous way such that all MPI ranks have the same random number. The argument is ignored unless it is of type MPIData, in which case a random number from the root rank is broadcasted to all MPI ranks. MPI syncronizing.

source
Rimu.DictVectors.walkernumberMethod
walkernumber(md::MPIData)

Compute the walkernumber of the distributed data on every MPI rank with MPI.Allreduce. MPI syncronizing.

source
Rimu.RMPI.freeMethod
free(obj::MPIOneSided)

De-reference the object, call finalizer and the garbage collector immediately. This is a syncronizing MPI call. Make sure that the object is not used later. MPI syncronizing.

source
Rimu.RMPI.mpi_combine_walkers!Method
mpi_combine_walkers!(target, source, [strategy])

Distribute the entries of source to the target data structure such that all entries in the target are on the process with the correct mpi rank as controlled by targetrank(). MPI syncronizing.

Note: the storage of the source is communicated rather than the source itself.

source
Rimu.RMPI.mpi_communicate_buffers!Method
mpi_communicate_buffers!(target::AbstractDVec{K,V}, buffers::Vector{<:Vector{V}})

Use MPI to communicate the contents of buffers and sort them into target. The length of buffers should be equal to mpi_size.

source
Rimu.RMPI.mpi_one_sidedFunction
mpi_one_sided(data, comm = mpi_comm(), root = mpi_root; capacity)

Declare data as mpi-distributed and set communication strategy to one-sided with remote memory access (RMA). capacity sets the capacity of the RMA windows.

Sets up the MPIData structure with MPIOneSided strategy.

source
Rimu.RMPI.mpi_seed_CRNGs!Function
mpi_seed_CRNGs!(seed = rand(Random.RandomDevice(), UInt))

Re-seed the random number generators in an MPI-safe way. If seed is provided, the random numbers from cRand() will follow a deterministic sequence.

Independence of the random number generators on different MPI ranks is achieved by adding hash(mpi_rank()) to seed.

source
Rimu.RMPI.putMethod
put(buf::Vector{T}, [len,] targetrank, s::MPIOneSided{T})

Deposit a vector buf into the MPI window s on rank targetrank. If len is given, only the first len elements are transmitted.

source
Rimu.RMPI.receive!Method
receive!(target, s::MPIPointToPoint, id)

Recieve from rank with id and move recieved values to target.

source
Rimu.RMPI.send!Method
send!(s::MPIPointToPoint{P})

Send the contents of the send buffers to all other ranks.

source
Rimu.RMPI.sort_and_count!Function
sort_and_count!(counts, displs, vec, order, (lo, hi), i_start, j_start)

Sort new spawns by target rank. While this is done, also count the values and calculate the offsets needed for MPI.Alltoallv!.

counts, displs, vec, and order are modified in-place. order should contain the values you want to sort vec by. (i.e. targetrank.(vec, s.np))

sort_and_count!(s::MPIAllToAll)

As above, but operating on the internal buffers of s. Note that s.targets is expected to contain the correct values to sort by.

source
Rimu.RMPI.@mpi_rootMacro
@mpi_root expr

Evaluate expression only on the root rank. Extra care needs to be taken as expr must not contain any code that involves syncronising MPI operations, i.e. actions that would require syncronous action of all MPI ranks.

Example:

wn = walkernumber(dv)   # an MPI syncronising function call that gathers
                        # information from all MPI ranks
@mpi_root @info "The current walker number is" wn # print info message on root only
source

Index