API
Rimu
Rimu.Rimu — ModuleRimuRandom Integrator for Many-Body Quantum Systems
Rimu.AllOverlaps — TypeAllOverlaps(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).
Rimu.ConstantTimeStep — TypeConstantTimeStep <: TimeStepStrategyKeep dτ constant.
Rimu.DeltaMemory — TypeDeltaMemory(Δ::Int) <: MemoryStrategyBefore 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̃>Rimu.DeltaMemory2 — TypeDeltaMemory2(Δ::Int) <: MemoryStrategyBefore 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.
Rimu.DeltaMemory3 — TypeDeltaMemory3(Δ::Int, level::Float64) <: MemoryStrategyBefore 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*rRimu.DontUpdate — TypeDontUpdate(; targetwalkers = 1_000_000) <: ShiftStrategyDon't update the shift. Return when targetwalkers is reached.
Rimu.DoubleLogProjected — TypeDoubleLogProjected(; target, projector, ζ = 0.08, ξ = ζ^2/4) <: ShiftStrategyStrategy 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)\]
Rimu.DoubleLogSumUpdate — TypeDoubleLogSumUpdate(; targetwalkers = 1000, ζ = 0.08, ξ = ζ^2/4, α = 1/2) <: ShiftStrategyStrategy 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/ζ.
Rimu.DoubleLogUpdate — TypeDoubleLogUpdate(; targetwalkers = 1000, ζ = 0.08, ξ = ζ^2/4) <: ShiftStrategyStrategy 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/ζ.
Rimu.DoubleLogUpdateAfterTargetWalkers — TypeDoubleLogUpdateAfterTargetWalkers(targetwalkers, ζ = 0.08, ξ = 0.0016) <: ShiftStrategyStrategy for updating the shift: After targetwalkers is reached, update the shift according to the log formula with damping parameter ζ and ξ. See DoubleLogUpdate.
Rimu.FciqmcRunStrategy — Type FciqmcRunStrategy{T}Abstract type representing the strategy for running and terminating fciqmc!(). The type parameter T is relevant for reporting the shift and the norm.
Implemented strategies:
Rimu.LogUpdate — TypeLogUpdate(ζ = 0.08) <: ShiftStrategyStrategy 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)\]
Rimu.LogUpdateAfterTargetWalkers — TypeLogUpdateAfterTargetWalkers(targetwalkers, ζ = 0.08) <: ShiftStrategyStrategy for updating the shift: After targetwalkers is reached, update the shift according to the log formula with damping parameter ζ. See LogUpdate.
Rimu.MemoryStrategy — TypeAbstract type for defining the strategy for injectimg memory noise. Implemented strategies:
Rimu.MultiScalar — TypeMultiScalarWrapper 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.
Rimu.NoMemory — TypeNoMemory <: MemoryStrategyDefault strategy for MemoryStrategy indicating that no memory noise will be used.
Rimu.NoStats — TypeNoStats(N=1) <: ReplicaStrategy{N}The default ReplicaStrategy. N replicas are run, but no statistics are collected.
Rimu.PostStepStrategy — TypePostStepStrategySubtypes of PostStepStrategy can be used to perform arbitrary computation on a replica after an FCIQMC step is finished and report the results.
Implemented strategies:
Note: A tuple of multiple strategies can be passed to lomc!. In that case, all reported column names must be distinct.
Interface:
A subtype of this type must implement post_step(::PostStepStrategy, ::ReplicaState).
Rimu.ProjectedEnergy — TypeProjectedEnergy(hamiltonian, projector; hproj=:vproj, vproj=:vproj) <: PostStepStrategyAfter 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.
Rimu.ProjectedMemory — TypeProjectedMemory(Δ::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.
Rimu.Projector — TypeProjector(name=projector) <: PostStepStrategyAfter each step, compute dot(projector, dv) and report it in the DataFrame under name. projector can be an AbstractDVec, or an AbstractProjector.
Rimu.PurgeNegatives — TypePurgeNegatives <: MemoryStrategyPurge all negative sign walkers.
Rimu.QMCState — TypeQMCStateHolds all information needed to run FCIQMC, except the data frame. Holds a NTuple of ReplicaStates and various strategies that control the algorithm.
Rimu.ReplicaState — TypeReplicaState(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 (seewalkernumber).params: theFCIQMCRunStrategy.id: appended to reported columns.
Rimu.ReplicaStrategy — TypeReplicaStrategy{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:
replica_stats(::ReplicaStrategy{N}, ::NTuple{N,ReplicaState})- return a tuple ofStrings orSymbolsof replica statistic names and a tuple of the values. These will be reported to theDataFramereturned bylomc!
Concrete implementations:
NoStats: run (possibly one) replica(s), but don't report any additional info.AllOverlaps: report overlaps between all pairs of replica vectors.
Rimu.Report — TypeReportInternal structure that holds the temporary reported values. See report!.
Rimu.ReportDFAndInfo — TypeReportDFAndInfo(; k=1, i=100, io=stdout, writeinfo=false) <: ReportingStrategyThe 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().
Rimu.ReportToFile — TypeReportToFile(; kwargs...) <: ReportingStrategyReporting 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. ADataFrameof this size is collected in memory and written to disk. When saving, an info message is also printed toio.save_if = true: if this value is true, save the report, otherwise ignore it. Usesave_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 emptyDataFrameis returned.io=stdout: TheIOto print messages to. Set todevnullif you don't want to see messages printed out.
Rimu.ReportingStrategy — TypeReportingStrategyAbstract type for strategies for reporting data in a DataFrame with report!().
Implemented strategies:
Interface:
A ReportingStrategy can define any of the following:
Rimu.RunTillLastStep — TypeRunTillLastStep(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
) <: FciqmcRunStrategyParameters for running fciqmc!() for a fixed number of time steps. For alternative strategies, see FciqmcRunStrategy.
Rimu.ShiftMemory — TypeShiftMemory(Δ::Int) <: MemoryStrategyEffectively replaces the fluctuating shift update procedure for the coefficient vector by an averaged shift over Δ timesteps, where Δ = 1 means no averaging.
Rimu.ShiftStrategy — TypeAbstract type for defining the strategy for updating the shift with update_shift(). Implemented strategies:
DontUpdateDoubleLogUpdate- default inlomc!()LogUpdateLogUpdateAfterTargetWalkers- FCIQMC standardDoubleLogUpdateAfterTargetWalkers
Rimu.SignCoherence — TypeSignCoherence(reference[; name=:coherence]) <: PostStepStrategyAfter 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.
Rimu.TimeStepStrategy — TypeTimeStepStrategyAbstract type for strategies for updating the time step with update_dτ(). Implemented strategies:
Rimu.TripleLogUpdate — TypeTripleLogUpdate(; targetwalkers = 1000, ζ = 0.08, ξ = ζ^2/4, η = 0.01) <: ShiftStrategyStrategy 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/ζ.
Rimu.WalkerLoneliness — TypeWalkerLoneliness(threshold=1) <: PostStepStrategyAfter each step, compute the proportion of configurations that are occupied by at most threshold walkers. Reports to a column named loneliness.
Rimu._n_walkers — Method_n_walkers(v, s_strat)Returns an estimate of the expected number of walkers as an integer.
Rimu.advance! — Methodadvance!(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.
Rimu.all_overlaps — Methodall_overlaps(operators, vectors)Get all overlaps between vectors and operators. This function is overlpaded for MPIData.
Rimu.apply_memory_noise! — Methodr = 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.
Rimu.fciqmc_col! — Methodfciqmc_col!(w, ham, add, num, shift, dτ)
fciqmc_col!(::Type{T}, args...)
-> spawns, deaths, clones, antiparticles, annihilationsSpawning 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:
IsDeterministic()deteministic algorithmIsStochasticInteger()stochastic version where the changes added toware purely integer, according to the FCIQMC algorithmIsStochasticWithThreshold(c)stochastic algorithm with floating point walkers.
Rimu.fciqmc_step! — Functionfciqmc_step!(Ĥ, v, shift, dτ, pnorm, w;
m_strat::MemoryStrategy = NoMemory()) -> ṽ, w̃, statsPerform 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 w̃, and the array stats = [spawns, deaths, clones, antiparticles, annihilations]. Stats will contain zeros when running in deterministic mode.
Rimu.finalize_report! — Methodfinalize_report!(::ReportingStrategy, report)Finalize the report. This function is called after all steps in lomc! have finished.
Rimu.lomc! — Methodlomc!(ham::AbstractHamiltonian, [v]; kwargs...) -> df, state
lomc!(state::QMCState, [df]; kwargs...) -> df, stateLinear 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, seeFciqmcRunStrategy; is mutatedlaststep- can be used to override information otherwise contained inparamss_strat::ShiftStrategy = DoubleLogUpdate(targetwalkers = 100, ζ = 0.08, ξ = ζ^2/4)- how to update theshift, seeShiftStrategymaxlength = 2 * s_strat.targetwalkers + 100- upper limit on the length ofv; when reached,lomc!will abortstyle = IsStochasticInteger()- setStochasticStylefor defaultv; unused ifvis specified.post_step::NTuple{N,<:PostStepStrategy} = ()- extract observables (e.g.ProjectedEnergy), seePostStepStrategy.replica::ReplicaStrategy = NoStats(1)- run several synchronised simulation, seeReplicaStrategy.r_strat::ReportingStrategy = ReportDFAndInfo()- how and when to report results, seeReportingStrategyτ_strat::TimeStepStrategy = ConstantTimeStep()- adjust time step dynamically, seeTimeStepStrategym_strat::MemoryStrategy = NoMemory()- experimental: inject memory noise, seeMemoryStrategythreading = :auto- can be used to control the use of multithreading (overridden bywm):auto- use multithreading ifs_strat.targetwalkers ≥ 500true- use multithreading if available (set shell variableJULIA_NUM_THREADS!)false- run on single thread
wm- working memory; if set, it controls the use of multithreading and overridesthreading; is mutateddf = DataFrame()- when called withAbstractHamiltonianargument, aDataFramecan be passed intolomc!that will be pushed into.
Return values
lomc! returns a named tuple with the following fields:
df: aDataFramewith all statistics being reported.state: aQMCStatethat 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)Rimu.post_step — Functionpost_step(::PostStepStrategy, ::ReplicaState) -> kvpairsCompute statistics after FCIQMC step. Should return a tuple of :key => value pairs. See also PostStepStrategy.
Rimu.refine_r_strat — Methodrefine_r_strat(r_strat::ReportingStrategy) -> r_stratInitialize the reporting strategy. This can be used to set up filenames or other attributes that need to be unique for a run of FCIQMC.
Rimu.replica_stats — Functionreplica_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.
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...).
Rimu.report! — Methodreport!(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.
Rimu.report_after_step — Methodreport_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.
Rimu.sort_into_targets! — Methodsort_into_targets!(target, source, stats) -> agg, wm, agg_statsAggregate 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.
Rimu.step_stats — Methodstep_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!.
Rimu.threshold_projected_deposit! — Methodthreshold_projected_deposit!This function performs threshold projection before spawning, but only for IsDynamicSemistochastic with the project_later parameter set to false.
Rimu.update_dvec! — Methodupdate_dvec!([::StochasticStyle,] dvec) -> dvec, ntPerform 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!
Rimu.update_dτ — Methodupdate_dτ(s<:TimeStepStrategy, dτ, tnorm) -> new dτUpdate the time step according to the strategy s.
Rimu.update_shift — Functionupdate_shift(s <: ShiftStrategy, shift, shiftMode, tnorm, pnorm, dτ, step, df, v_new, v_old)Update the shift according to strategy s. See ShiftStrategy.
Reexported Submodules
Hamiltonians
Link to Module Hamiltionians.jl
BitStringAddresses
Link to Module BitStringAddresses.jl
DictVectors
Rimu.DictVectors — ModuleModule 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
Rimu.DictVectors.AbstractDVec — TypeDictVectors.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:
StochasticStylestorage(dv)returns anAbstractDictstoring the raw data with possibly differentvaltypethanV.
Rimu.DictVectors.AbstractProjector — TypeAbstract supertype for projectors to be used in in lieu of DVecs or Vectors.
Rimu.DictVectors.CoherentInitiator — TypeCoherentInitiator(threshold) <: InitiatorRuleInitiator 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.
Rimu.DictVectors.DVec — TypeDVec{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 aDVecwithdictfor storage. Note that the data may or may not be copied.DVec(args...[; style, capacity]):args...are passed to theDictconstructor. TheDictis used for storage.DVec{K,V}([; style, capacity]): create an emptyDVec{K,V}.DVec(dv::AbstractDVec[; style, capacity]): create aDVecwith the same contents asadv. Thestyleis inherited fromdvby 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.0Rimu.DictVectors.FrozenDVec — TypeFrozenDVecSee: freeze.
Rimu.DictVectors.Initiator — TypeInitiator(threshold) <: InitiatorRuleInitiator 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.
Rimu.DictVectors.InitiatorDVec — TypeInitiatorDVec{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 anInitiatorDVecwithdictfor storage. Note that the data may or may not be copied.InitiatorDVec(args...[; style, initiator, capacity]):args...are passed to theDictconstructor. TheDictis used for storage.InitiatorDVec{K,V}([; style, initiator, capacity]): create an emptyInitiatorDVec{K,V}.InitiatorDVec(dv::AbstractDVec[; style, initiator, capacity]): create anInitiatorDVecwith the same contents asdv. Thestyleis inherited fromdvby default.
Keyword arguments
style: A validStochasticStyle. The default is selected based on theInitiatorDVec'svaltype(seedefault_style). If a style is given and thevaltypedoes not match thestyle'seltype, the values are converted to an appropriate type.initiator = Initiator(1): A validInitiatorRule. SeeInitiator.capacity: Indicative size asInt. Optional. Sets the initial size of theInitiatorDVecviasizehint!.
Rimu.DictVectors.InitiatorIterator — TypeInitiatorIteratorIterator over pairs or values of an InitiatorDVec. Supports the SplittablesBase interface.
Rimu.DictVectors.InitiatorRule — TypeInitiatorRule{V}Abstract type for defining initiator rules for InitiatorDVec. Concrete implementations:
When defining a new InitiatorRule, also define a corresponding method for value!
Rimu.DictVectors.InitiatorValue — TypeInitiatorValue{V}(; safe::V, unsafe::V, initiator::V) where VComposite "walker" with three fields. For use with InitiatorDVecs.
Rimu.DictVectors.IsDeterministic — TypeIsDeterministic{T=Float64}() <: StochasticStyle{T}Trait for generalised vector of configuration indicating deterministic propagation of walkers.
See also StochasticStyle.
Rimu.DictVectors.IsDynamicSemistochastic — TypeIsDynamicSemistochastic{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.gabs_threshold = 0.1 * target_walkers.proj_threshold = 1.0: Values below this number are stochastically projected to this value or zero. See alsoIsStochasticWithThreshold.
See also StochasticStyle.
Rimu.DictVectors.IsStochastic2Pop — TypeIsStochastic2Pop{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.
Rimu.DictVectors.IsStochasticInteger — TypeIsStochasticInteger{T=Int}() <: StochasticStyle{T}Trait for generalised vector of configurations indicating stochastic propagation as seen in the original FCIQMC algorithm.
See also StochasticStyle.
Rimu.DictVectors.IsStochasticWithThreshold — TypeIsStochasticWithThreshold(threshold=1.0) <: StochasticStyleTrait 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.
Rimu.DictVectors.Norm1ProjectorPPop — TypeNorm1ProjectorPPop()Results in computing the one-norm per population when used in dot(). E.g.
dot(Norm1ProjectorPPop(),x)
-> norm(real.(x),1) + im*norm(imag.(x),1)See also ReportingStrategy for use of projectors in FCIQMC.
Rimu.DictVectors.Norm2Projector — TypeNorm2Projector()Results in computing the two-norm when used in dot(). E.g.
dot(NormProjector(),x)
-> norm(x,2) # with type Float64See also ReportingStrategy for use of projectors in FCIQMC.
Rimu.DictVectors.NormProjector — TypeNormProjector()Results in computing the one-norm when used in dot(). E.g.
dot(NormProjector(),x)
-> norm(x,1)NormProjector() thus represents the vector sign.(x).
See also ReportingStrategy for use of projectors in FCIQMC.
Rimu.DictVectors.PopsProjector — TypePopsProjector()Results in computing the projection of one population on the other when used in dot(). E.g.
dot(PopsProjector(),x)
-> real(x) ⋅ imag(x)See also ReportingStrategy for use of projectors in FCIQMC.
Rimu.DictVectors.SimpleInitiator — TypeSimpleInitiator(threshold) <: InitiatorRuleSimplified 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.
Rimu.DictVectors.StochasticStyle — TypeStochasticStyle(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
IsStochasticInteger- integer walker FCIQMCIsDeterministic- perform deterministic variant of power methodIsStochasticWithThreshold- floating point walker FCIQMCIsDynamicSemistochastic
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.
Rimu.DictVectors.StyleUnknown — TypeStyleUnknown{T}() <: StochasticStyleTrait for value types not (currently) compatible with FCIQMC. This style makes it possible to construct dict vectors with unsupported valtypes.
See also StochasticStyle.
Rimu.DictVectors.UniformProjector — TypeUniformProjector()Represents a vector with all elements 1. To be used with dot(). Minimizes memory allocations.
UniformProjector()⋅v == sum(v)
dot(UniformProjector(), LO, v) == sum(LO*v)See also ReportingStrategy for use of projectors in FCIQMC.
OrderedCollections.freeze — Methodfreeze(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.
Rimu.DictVectors.add! — Methodadd!(x::AbstractDVec,y::AbstactDVec)Inplace add x+y and store result in x.
Rimu.DictVectors.default_style — Methoddefault_style(::Type)Pick a StochasticStyle based on the value type. Throws an error if no known default style is known.
Rimu.DictVectors.deposit! — Methoddeposit!(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.
Rimu.DictVectors.deposit! — Methoddeposit!(w::InitiatorDVec, add, val, p_add=>p_val)Add val into w at address add as an InitiatorValue.
Rimu.DictVectors.localpart — Methodlocalpart(dv) -> AbstractDVecGet the part of dv that is located on this MPI rank. Returns dv itself for DictVectors.
Rimu.DictVectors.storage — Functionstorage(dvec) -> AbstractDictReturn the raw storage associated with dvec as an AbstractDict. Used in MPI communication.
Rimu.DictVectors.value — Functionvalue(i::InitiatorRule, v::InitiatorValue)Convert the InitiatorValue v into a scalar value according to the InitiatorRule i.
Internal function that implements functionality of InitiatorDVec.
Rimu.DictVectors.walkernumber — Methodwalkernumber(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.
Rimu.DictVectors.zero! — Methodzero!(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.
ConsistentRNG
Link to Module ConsistentRNG.jl
Unexported Submodules
StatsTools
Link to Module Rimu/StatsTools
Blocking
Rimu.Blocking — ModuleBlocking
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)).
Rimu.Blocking.autoblock — Methodautoblock(df::DataFrame; start = 1, stop = size(df)[1])
-> s̄, σs, ē, σe, kDetermine mean shift s̄ 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.
Rimu.Blocking.autoblock — Methodautoblock(dftup::Tuple; start = 1, stop = size(dftup[1])[1])
-> s̄1, σs1, s̄2, σs2, ē1, σe1, ē2, σe2, ēH, σeH, kReplica 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.
Rimu.Blocking.autocovariance — Methodautocovariance(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).
Rimu.Blocking.blockAndMTest — Methodv̄, σ, σσ, k, df = blockAndMTest(v::Vector)Perform a blocking analysis and M-test on v returning the mean v̄, standard error σ, its error σσ, the number of blocking steps k, and the DataFrame df with blocking data.
Rimu.Blocking.blocker — Methodblocker(v::Vector) -> new_v::VectorReblock 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.
Rimu.Blocking.blocking — Methodblocking(x::Vector,y::Vector) -> df::DataFramePerform 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 datasetx;mean_y,SD_y,SE_y,SE_SE_y= ditto. for datasety;Covariance= the covariance between data inxandy;mean_f=x̄/ȳ;SE_f= standard error estimated forx̄/ȳ.
Rimu.Blocking.blocking — Methodblocking(v::Vector; corrected::Bool=true) -> dfPerform 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).
Rimu.Blocking.blocking_old — Methodblocking(v::Vector; typos = nothing) -> dfPerform 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- calculateMand standard error as written in Jonsson.
Rimu.Blocking.combination_division — Methodcombination_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).
Rimu.Blocking.covariance — Methodcovariance(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).
Rimu.Blocking.gW — MethodgW(norm::AbstractArray, shift::AbstractArray, dt, [b]; pad = :true) -> g
gW(df::DataFrame, [b]; pad = :true) -> gCompute 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.
Rimu.Blocking.growthWitness — MethodgrowthWitness(norm::AbstractArray, shift::AbstractArray, dt; b = 30, pad = :true) -> g
growthWitness(df::DataFrame; b = 30, pad = :true) -> gCompute 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 S̄ 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.
Rimu.Blocking.mtest — Methodmtest(df::DataFrame; warn = true) -> kThe "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.
Rimu.Blocking.se — Methodse(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).
Rimu.Blocking.smoothen — Methodsmoothen(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.
RMPI
Rimu.RMPI — ModuleModule for providing MPI functionality for Rimu.
Rimu.RMPI.mpi_root — ConstantDefault MPI root for RMPI.
Rimu.RMPI.MPIAllToAll — Type MPIAllToAllAll-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 typePonnpprocesses with current rankid.
Rimu.RMPI.MPIData — TypeMPIData(data; kwargs...)Wrapper used for signaling that this data is part of a distributed data structure and communication should happen with MPI.
Keyword arguments:
setup = mpi_point_to_point- controls the communication stratgympi_one_sideduses one-sided communication with remote memory access (RMA), setsMPIOneSidedstrategy.mpi_point_to_pointusesMPIPointTOPointstrategy.mpi_all_to_allusesMPIAllToAllstrategy.mpi_no_exchangesetsMPINoWalkerExchangestrategy. Experimental. Use with caution!
comm = mpi_comm()root = mpi_root- The rest of the keyword arguments are passed to
setup.
Rimu.RMPI.MPIDataIterator — TypeMPIDataIterator{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.
Rimu.RMPI.MPINoWalkerExchange — TypeMPINoWalkerExchange(nprocs, my_rank, comm)Strategy for for not exchanging walkers between ranks. Consequently there will be no cross-rank annihilations.
Rimu.RMPI.MPIOneSided — TypeMPIOneSided(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().
Rimu.RMPI.MPIPointToPoint — TypeMPIPointToPoint{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 typePonnpprocesses with current rankid.
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.
Base.length — Methodlength(md::MPIData)Compute the length of the distributed data on every MPI rank with MPI.Allreduce. MPI syncronizing.
LinearAlgebra.norm — Functionnorm(md::MPIData, p=2)Compute the norm of the distributed data on every MPI rank with MPI.Allreduce. MPI syncronizing.
Rimu.ConsistentRNG.check_crng_independence — MethodConsistentRNGs.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.
Rimu.ConsistentRNG.sync_cRandn — Methodsync_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.
Rimu.DictVectors.walkernumber — Methodwalkernumber(md::MPIData)Compute the walkernumber of the distributed data on every MPI rank with MPI.Allreduce. MPI syncronizing.
Rimu.RMPI.copy_to_local! — Methodcopy_to_local!(target, md::MPIData)Collect all pairs in md from all ranks and store them in target. In-place version of copy_to_local.
Rimu.RMPI.copy_to_local — Methodcopy_to_local(md::MPIData)Collect all pairs in md from all ranks and store them in a local AbstractDVec.
Rimu.RMPI.free — Methodfree(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.
Rimu.RMPI.is_mpi_root — Functionis_mpi_root(root = mpi_root)Returns true if called from the root rank
Rimu.RMPI.mpi_all_to_all — Functionmpi_all_to_all(data, comm = mpi_comm(), root = mpi_root)Declare data as mpi-distributed and set communication strategy to all-to-all.
Sets up the MPIData structure with MPIAllToAll strategy.
Rimu.RMPI.mpi_allprintln — Methodmpi_allprintln(args...)Print a message to stdout from each rank separately, in order. MPI synchronizing.
Rimu.RMPI.mpi_barrier — Functionmpi_barrier(comm = mpi_comm())The MPI barrier with optional argument. MPI syncronizing.
Rimu.RMPI.mpi_combine_walkers! — Methodmpi_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.
Rimu.RMPI.mpi_comm — MethodDefault MPI communicator for RMPI.
Rimu.RMPI.mpi_communicate_buffers! — Methodmpi_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.
Rimu.RMPI.mpi_no_exchange — Functionmpi_no_exchange(data, comm = mpi_comm(), root = mpi_root)Declare data as mpi-distributed and set communication strategy to MPINoWalkerExchange. Sets up the MPIData structure with MPINoWalkerExchange strategy.
Rimu.RMPI.mpi_one_sided — Functionmpi_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.
Rimu.RMPI.mpi_point_to_point — Functionmpi_point_to_point(data, comm = mpi_comm(), root = mpi_root)Declare data as mpi-distributed and set communication strategy to point-to-point.
Sets up the MPIData structure with MPIPointToPoint strategy.
Rimu.RMPI.mpi_rank — Functionmpi_rank(comm = mpi_comm())Return the current MPI rank.
Rimu.RMPI.mpi_seed_CRNGs! — Functionmpi_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.
Rimu.RMPI.mpi_size — Functionmpi_size(comm = mpi_comm())Size of MPI communicator.
Rimu.RMPI.mpi_synchronize! — Methodmpi_synchronize!(md::MPIData)Synchronize md, ensuring its contents are distributed among ranks correctly.
Rimu.RMPI.next_mpiID — Functionnext_mpiID()Produce a new ID number for MPI distributed objects. Uses an internal counter.
Rimu.RMPI.put — Methodput(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.
Rimu.RMPI.receive! — Methodreceive!(target, s::MPIPointToPoint, id)Recieve from rank with id and move recieved values to target.
Rimu.RMPI.recvbuff — Methodrecvbuff(s::MPIPointToPoint)Get the receive buffer.
Rimu.RMPI.send! — Methodsend!(s::MPIPointToPoint{P})Send the contents of the send buffers to all other ranks.
Rimu.RMPI.sendbuff — Methodsendbuff(s::MPIPointToPoint, id)Get the send buffer associated with id.
Rimu.RMPI.sort_and_count! — Functionsort_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.
Rimu.RMPI.targetrank — Methodtargetrank(key, np)Compute the rank where the key belongs.
Rimu.RMPI.@mpi_root — Macro@mpi_root exprEvaluate 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 onlyRimu.RMPI.@swap! — Macro@swap! arr i jSwap the i-th and j-th indices in arr.
Index
Rimu.BlockingRimu.ConsistentRNGRimu.DictVectorsRimu.HamiltoniansRimu.RMPIRimu.RimuRimu.StatsToolsRimu.ConsistentRNG.CRNGsRimu.RMPI.mpi_rootMonteCarloMeasurements.ParticlesRimu.AllOverlapsRimu.BitStringAddresses.AbstractFockAddressRimu.BitStringAddresses.BitStringRimu.BitStringAddresses.BoseFSRimu.BitStringAddresses.BoseFS2CRimu.BitStringAddresses.CompositeFSRimu.BitStringAddresses.FermiFSRimu.BitStringAddresses.SingleComponentFockAddressRimu.ConsistentRNG.CRNGRimu.ConstantTimeStepRimu.DeltaMemoryRimu.DeltaMemory2Rimu.DeltaMemory3Rimu.DictVectors.AbstractDVecRimu.DictVectors.AbstractProjectorRimu.DictVectors.CoherentInitiatorRimu.DictVectors.DVecRimu.DictVectors.FrozenDVecRimu.DictVectors.InitiatorRimu.DictVectors.InitiatorDVecRimu.DictVectors.InitiatorIteratorRimu.DictVectors.InitiatorRuleRimu.DictVectors.InitiatorValueRimu.DictVectors.IsDeterministicRimu.DictVectors.IsDynamicSemistochasticRimu.DictVectors.IsStochastic2PopRimu.DictVectors.IsStochasticIntegerRimu.DictVectors.IsStochasticWithThresholdRimu.DictVectors.Norm1ProjectorPPopRimu.DictVectors.Norm2ProjectorRimu.DictVectors.NormProjectorRimu.DictVectors.PopsProjectorRimu.DictVectors.SimpleInitiatorRimu.DictVectors.StochasticStyleRimu.DictVectors.StyleUnknownRimu.DictVectors.UniformProjectorRimu.DontUpdateRimu.DoubleLogProjectedRimu.DoubleLogSumUpdateRimu.DoubleLogUpdateRimu.DoubleLogUpdateAfterTargetWalkersRimu.FciqmcRunStrategyRimu.Hamiltonians.AbstractHamiltonianRimu.Hamiltonians.BoseHubbardMom1D2CRimu.Hamiltonians.BoseHubbardReal1D2CRimu.Hamiltonians.ExtendedHubbardReal1DRimu.Hamiltonians.HardwallBoundariesRimu.Hamiltonians.HubbardMom1DRimu.Hamiltonians.HubbardReal1DRimu.Hamiltonians.HubbardRealSpaceRimu.Hamiltonians.LOStructureRimu.Hamiltonians.LadderBoundariesRimu.Hamiltonians.LatticeGeometryRimu.Hamiltonians.MatrixHamiltonianRimu.Hamiltonians.PeriodicBoundariesRimu.LogUpdateRimu.LogUpdateAfterTargetWalkersRimu.MemoryStrategyRimu.MultiScalarRimu.NoMemoryRimu.NoStatsRimu.PostStepStrategyRimu.ProjectedEnergyRimu.ProjectedMemoryRimu.ProjectorRimu.PurgeNegativesRimu.QMCStateRimu.RMPI.MPIAllToAllRimu.RMPI.MPIDataRimu.RMPI.MPIDataIteratorRimu.RMPI.MPINoWalkerExchangeRimu.RMPI.MPIOneSidedRimu.RMPI.MPIPointToPointRimu.ReplicaStateRimu.ReplicaStrategyRimu.ReportRimu.ReportDFAndInfoRimu.ReportToFileRimu.ReportingStrategyRimu.RunTillLastStepRimu.ShiftMemoryRimu.ShiftStrategyRimu.SignCoherenceRimu.StatsTools.BlockingResultRimu.StatsTools.RatioBlockingResultRimu.TimeStepStrategyRimu.TripleLogUpdateRimu.WalkerLonelinessBase.:*Base.lengthLinearAlgebra.normMeasurements.measurementOrderedCollections.freezeRimu.BitStringAddresses.chunk_bitsRimu.BitStringAddresses.chunk_typeRimu.BitStringAddresses.find_modeRimu.BitStringAddresses.find_occupied_modeRimu.BitStringAddresses.is_occupiedRimu.BitStringAddresses.move_particleRimu.BitStringAddresses.near_uniformRimu.BitStringAddresses.num_bitsRimu.BitStringAddresses.num_chunksRimu.BitStringAddresses.num_componentsRimu.BitStringAddresses.num_modesRimu.BitStringAddresses.num_occupied_modesRimu.BitStringAddresses.num_particlesRimu.BitStringAddresses.occupied_modesRimu.BitStringAddresses.onrRimu.BitStringAddresses.top_chunk_bitsRimu.Blocking.autoblockRimu.Blocking.autoblockRimu.Blocking.autocovarianceRimu.Blocking.blockAndMTestRimu.Blocking.blockerRimu.Blocking.blockingRimu.Blocking.blockingRimu.Blocking.blocking_oldRimu.Blocking.combination_divisionRimu.Blocking.covarianceRimu.Blocking.gWRimu.Blocking.growthWitnessRimu.Blocking.mtestRimu.Blocking.seRimu.Blocking.smoothenRimu.ConsistentRNG.cRandRimu.ConsistentRNG.cRandnRimu.ConsistentRNG.check_crng_independenceRimu.ConsistentRNG.check_crng_independenceRimu.ConsistentRNG.newChildRNGRimu.ConsistentRNG.seedCRNG!Rimu.ConsistentRNG.sync_cRandnRimu.ConsistentRNG.sync_cRandnRimu.ConsistentRNG.trngRimu.DictVectors.add!Rimu.DictVectors.default_styleRimu.DictVectors.deposit!Rimu.DictVectors.deposit!Rimu.DictVectors.localpartRimu.DictVectors.storageRimu.DictVectors.valueRimu.DictVectors.walkernumberRimu.DictVectors.walkernumberRimu.DictVectors.zero!Rimu.Hamiltonians.diagonal_elementRimu.Hamiltonians.dimensionRimu.Hamiltonians.get_offdiagonalRimu.Hamiltonians.neighbour_siteRimu.Hamiltonians.num_neighboursRimu.Hamiltonians.num_offdiagonalsRimu.Hamiltonians.offdiagonalsRimu.Hamiltonians.random_offdiagonalRimu.Hamiltonians.starting_addressRimu.RMPI.copy_to_localRimu.RMPI.copy_to_local!Rimu.RMPI.freeRimu.RMPI.is_mpi_rootRimu.RMPI.mpi_all_to_allRimu.RMPI.mpi_allprintlnRimu.RMPI.mpi_barrierRimu.RMPI.mpi_combine_walkers!Rimu.RMPI.mpi_commRimu.RMPI.mpi_communicate_buffers!Rimu.RMPI.mpi_no_exchangeRimu.RMPI.mpi_one_sidedRimu.RMPI.mpi_point_to_pointRimu.RMPI.mpi_rankRimu.RMPI.mpi_seed_CRNGs!Rimu.RMPI.mpi_sizeRimu.RMPI.mpi_synchronize!Rimu.RMPI.next_mpiIDRimu.RMPI.putRimu.RMPI.receive!Rimu.RMPI.recvbuffRimu.RMPI.send!Rimu.RMPI.sendbuffRimu.RMPI.sort_and_count!Rimu.RMPI.targetrankRimu.StatsTools.autocovarianceRimu.StatsTools.blockerRimu.StatsTools.blocking_analysisRimu.StatsTools.blocks_with_mRimu.StatsTools.growth_estimatorRimu.StatsTools.growth_witnessRimu.StatsTools.mean_and_seRimu.StatsTools.med_and_errsRimu.StatsTools.mixed_estimatorRimu.StatsTools.mtestRimu.StatsTools.pseudo_covRimu.StatsTools.ratio_estimatorsRimu.StatsTools.ratio_of_meansRimu.StatsTools.ratio_with_errsRimu.StatsTools.replica_fidelityRimu.StatsTools.smoothenRimu.StatsTools.to_measurementRimu.StatsTools.w_expRimu.StatsTools.w_linRimu.StatsTools.x_by_y_linearRimu._n_walkersRimu.advance!Rimu.all_overlapsRimu.apply_memory_noise!Rimu.fciqmc_col!Rimu.fciqmc_step!Rimu.finalize_report!Rimu.lomc!Rimu.post_stepRimu.refine_r_stratRimu.replica_statsRimu.report!Rimu.report!Rimu.report_after_stepRimu.sort_into_targets!Rimu.step_statsRimu.threshold_projected_deposit!Rimu.update_dvec!Rimu.update_dτRimu.update_shiftStatistics.covRimu.RMPI.@mpi_rootRimu.RMPI.@swap!