API
Rimu
Rimu.Rimu
— ModuleRimu
Random integrators for many-body quantum systems
Welcome to Rimu
version 0.10.3. Read the documentation online.
Rimu.PACKAGE_VERSION
— ConstantRimu.PACKAGE_VERSION
Constant that contains the current VersionNumber
of Rimu
.
Rimu.AllOverlaps
— TypeAllOverlaps(num_replicas=2; operator=nothing, transform=nothing, vecnorm=true) <: ReplicaStrategy{num_replicas}
Run num_replicas
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}_dot_c{j}
for vector-vector overlaps, and c{i}_Op{k}_c{j}
for operator overlaps.
See lomc!
, ReplicaStrategy
and AbstractHamiltonian
(for an interface for implementing operators).
Transformed Hamiltonians
If a transformed Hamiltonian G
has been passed to lomc!
then overlaps can be calculated by passing the same transformed Hamiltonian to AllOverlaps
by setting transform=G
. A warning is given if these two Hamiltonians do not match.
Implemented transformations are:
In the case of a transformed Hamiltonian the overlaps are defined as follows. For a similarity transformation G
of the Hamiltonian (see e.g. GutzwillerSampling
.)
\[ \hat{G} = f \hat{H} f^{-1}.\]
The expectation value of an operator $\hat{A}$ is
\[ \langle \hat{A} \rangle = \langle \psi | \hat{A} | \psi \rangle = \frac{\langle \phi | f^{-1} \hat{A} f^{-1} | \phi \rangle}{\langle \phi | f^{-2} | \phi \rangle}\]
where
\[ | \phi \rangle = f | \psi \rangle\]
is the (right) eigenvector of $\hat{G}$ and $| \psi \rangle$ is an eigenvector of $\hat{H}$.
For a K-tuple of input operators (\hat{A}_1, ..., \hat{A}_K)
, overlaps of $\langle \phi | f^{-1} \hat{A} f^{-1} | \phi \rangle$ are reported as c{i}_Op{k}_c{j}
. The correct vector-vector overlap $\langle \phi | f^{-2} | \phi \rangle$ is reported last as c{i}_Op{K+1}_c{j}
. This is in addition to the bare vector-vector overlap $\langle \phi | f^{-2} | \phi \rangle$ that is reported as c{i}_dot_c{j}
.
In either case the c{i}_dot_c{j}
overlap can be omitted with the flag vecnorm=false
.
Rimu.ConstantTimeStep
— TypeConstantTimeStep <: TimeStepStrategy
Keep dτ
constant.
Rimu.DontUpdate
— TypeDontUpdate(; targetwalkers = 1_000_000) <: ShiftStrategy
Don't update the shift
. Return when targetwalkers
is reached.
See ShiftStrategy
, lomc!
.
Rimu.DoubleLogProjected
— TypeDoubleLogProjected(; 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)\]
Note that adjusting the keyword maxlength
in lomc!
is advised as the default may not be appropriate.
See ShiftStrategy
, lomc!
.
Rimu.DoubleLogSumUpdate
— TypeDoubleLogSumUpdate(; 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/ζ.
See ShiftStrategy
, lomc!
.
Rimu.DoubleLogUpdate
— TypeDoubleLogUpdate(; 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/ζ.
See ShiftStrategy
, lomc!
.
Rimu.DoubleLogUpdateAfterTargetWalkers
— TypeDoubleLogUpdateAfterTargetWalkers(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
, ShiftStrategy
, lomc!
.
Rimu.FciqmcRunStrategy
— Type FciqmcRunStrategy{T}
Abstract type representing the strategy for running and terminating lomc!()
. The type parameter T
is relevant for reporting the shift and the norm.
Implemented strategies:
Rimu.FirstOrderTransitionOperator
— TypeFirstOrderTransitionOperator(hamiltonian, shift, dτ) <: AbstractHamiltonian
First order transition operator
\[𝐓 = 1 + dτ(S - 𝐇)\]
where $𝐇$ is the hamiltonian
and $S$ is the shift
.
$𝐓$ represents the first order expansion of the exponential evolution operator of the imaginary-time Schrödinger equation (Euler step) and repeated application will project out the ground state eigenvector of the hamiltonian
. It is the transition operator used in FCIQMC.
Rimu.LogUpdate
— TypeLogUpdate(ζ = 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)\]
See ShiftStrategy
, lomc!
.
Rimu.LogUpdateAfterTargetWalkers
— TypeLogUpdateAfterTargetWalkers(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
, ShiftStrategy
, lomc!
.
Rimu.MultiScalar
— TypeMultiScalar
Wrapper over a tuple that supports +
, *
, min
, and max
. Used with MPI communication because SVector
s 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 MultiScalar
s with the same types can be operated on. This is a feature, as it forces type stability.
Rimu.NoStats
— TypeNoStats(N=1) <: ReplicaStrategy{N}
The default ReplicaStrategy
. N
replicas are run, but no statistics are collected.
See also lomc!
.
Rimu.PostStepStrategy
— TypePostStepStrategy
Subtypes 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=:hproj, 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 projected_energy
. 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.
See also projected_energy
, ratio_of_means
, mixed_estimator
, and PostStepStrategy
.
Rimu.Projector
— TypeProjector(name=projector) <: PostStepStrategy
After each step, compute dot(projector, dv)
and report it in the DataFrame
under name
. projector
can be an AbstractDVec
, or an AbstractProjector
.
Rimu.QMCState
— TypeQMCState
Holds all information needed to run lomc!
, except the dataframe. Holds an NTuple
of ReplicaState
s, the Hamiltonian, and various strategies that control the algorithm. Constructed and returned by lomc!
.
Rimu.ReplicaState
— TypeReplicaState(v, wm, 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.pv
: vector from the previous step.wm
: working memory.pnorm
: previous walker number (seewalkernumber
).params
: theFciqmcRunStrategy
.id
: string ID appended to reported column names.
See also QMCState
, ReplicaStrategy
, replica_stats
, lomc!
.
Rimu.ReplicaStrategy
— TypeReplicaStrategy{N}
Supertype for strategies that can be passed to lomc!
and control how many replicas are used, and what information is computed and returned. The number of replicas is N
.
Concrete implementations
NoStats
: run (possibly one) replica(s), but don't report any additional info.AllOverlaps
: report overlaps between all pairs of replica vectors.
Interface
A subtype of ReplicaStrategy{N}
must implement the following function:
Rimu.replica_stats
- return a tuple ofString
s orSymbols
of names for replica statistics and a tuple of the values. These will be reported to theDataFrame
returned bylomc!
.
Rimu.Report
— TypeReport()
Internal structure that holds the temporary reported values as well as metadata.
See report!
, report_metadata!
, get_metadata
.
Rimu.ReportDFAndInfo
— TypeReportDFAndInfo(; reporting_interval=1, info_interval=100, io=stdout, writeinfo=false) <: ReportingStrategy
The default ReportingStrategy
. Report every reporting_interval
th step to a DataFrame
and write info message to io
every info_interval
th reported 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...) <: ReportingStrategy
ReportingStrategy
that writes the report directly to a file in the Arrow
format. Useful when dealing with long jobs or large numbers of replicas, when the report can incur a significant memory cost.
The arrow file can be read back in with load_df(filename)
or using Arrow; Arrow.Table(filename)
.
Keyword arguments
filename = "out.arrow"
: the file to report to. If the file already exists, a new file is created.reporting_interval = 1
: interval between simulation steps that are reported.chunk_size = 1000
: the size of each chunk that is written to the file. ADataFrame
of this size is collected in memory and written to disk. When saving, an info message is also printed toio
.save_if =
is_mpi_root()
: if this value is true, save the report, otherwise ignore it.return_df = false
: if this value is true, read the file and return the data frame at the end of computation. Otherwise, an emptyDataFrame
is returned.io = stdout
: TheIO
to print messages to. Set todevnull
if you don't want to see messages printed out.compress = :zstd
: compression algorithm to use. Can be:zstd
,:lz4
ornothing
.
Rimu.ReportingStrategy
— TypeReportingStrategy
Abstract 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
) <: FciqmcRunStrategy
Parameters for running lomc!()
for a fixed number of time steps. For alternative strategies, see FciqmcRunStrategy
.
Rimu.ShiftStrategy
— TypeAbstract type for defining the strategy for updating the shift
. Passed as a parameter to lomc!
.
Implemented strategies:
DontUpdate
DoubleLogUpdate
- default inlomc!()
LogUpdate
LogUpdateAfterTargetWalkers
- FCIQMC standardDoubleLogUpdateAfterTargetWalkers
Rimu.SignCoherence
— TypeSignCoherence(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
.
Rimu.SingleParticleDensity
— TypeSingleParticleDensity(; save_every=1, component) <: PostStepStrategy
PostStepStrategy
to compute the diagonal single_particle_density
. It records a Tuple
with the same eltype
as the vector.
Computing the density at every time step can be expensive. This cost can be reduced by setting the save_every
argument to a higher value. If the value is set, a vector of zeros is recorded when the saving is skipped.
If the address type has multiple components, the component
argument can be used to compute the density on a per-component basis.
The density is not normalized, and must be divided by the vector norm(⋅,2)
squared.
See also
Rimu.TimeStepStrategy
— TypeTimeStepStrategy
Abstract type for strategies for updating the time step with update_dτ()
. Implemented strategies:
Rimu.Timer
— TypeTimer <: PostStepStrategy
Record current time after every step. See Base.Libc.time
for information on what time is recorded.
Rimu.TripleLogUpdate
— TypeTripleLogUpdate(; 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/ζ.
See ShiftStrategy
, lomc!
.
Rimu.WalkerLoneliness
— TypeWalkerLoneliness(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
.
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, working_memories, vecnorm=true)
Get all overlaps between vectors and operators. This function is overloaded for MPIData
. The flag vecnorm
can disable the vector-vector overlap c{i}_dot_c{j}
.
Rimu.check_transform
— Methodcheck_transform(r::AllOverlaps, ham)
Check that the transformation provided to r::AllOverlaps
matches the given Hamiltonian ham
. Used as a sanity check before starting main lomc!
loop.
Rimu.default_logger
— Methoddefault_logger(args...)
Reset the global_logger
to Logging.ConsoleLogger
. Undoes the effect of smart_logger
. Arguments are passed on to Logging.ConsoleLogger
.
Rimu.default_starting_vector
— Methoddefault_starting_vector(hamiltonian::AbstractHamiltonian; kwargs...)
default_starting_vector(
address=starting_address(hamiltonian);
style=IsStochasticInteger(),
initiator=NonInitiator(),
threading=nothing
)
Return a default starting vector for lomc!
. The default choice for the starting vector is
v = PDVec(address => 10; style, initiator)
if threading is available, or otherwise
v = DVec(address => 10; style)
if initiator == NonInitiator()
, and
v = InitiatorDVec(address => 10; style, initiator)
if not. See PDVec
, DVec
, InitiatorDVec
, StochasticStyle
, and [InitiatorRule
].
Rimu.finalize_report!
— Methodfinalize_report!(::ReportingStrategy, report)
Finalize the report. This function is called after all steps in lomc!
have finished.
Rimu.get_metadata
— Methodget_metadata(report::Report, key)
Get metadata key
from report
. key
is converted to a String
.
See also report_metadata!
, Report
, report!
.
Rimu.lomc!
— Methodlomc!(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
. The details of the simulation are controlled by the optional keyword arguments and by the type of the optional starting vector v
. Alternatively, a QMCState
can be passed in to continue a previous simulation.
Common keyword arguments and defaults:
laststep = 100
- controls the number of steps.dτ = 0.01
- time step.targetwalkers = 1000
- target for the 1-norm of the coefficient vector.address = starting_address(ham)
- set starting address for defaultv
andshift
.style = IsStochasticInteger()
- setStochasticStyle
for defaultv
; unused ifv
is specified.initiator = NonInitiator()
- setInitiatorRule
for defaultv
; unused ifv
is specified.threading
- default is to use multithreading and MPI if multiple threads are available. Set totrue
to forcePDVec
for the starting vector,false
for serial computation; unused ifv
is specified.shift = diagonal_element(ham, address)
- initial value of shift.post_step::NTuple{N,<:PostStepStrategy} = ()
- extract observables (e.g.ProjectedEnergy
), seePostStepStrategy
.replica::ReplicaStrategy = NoStats(1)
- run several synchronised simulations, seeReplicaStrategy
.r_strat::ReportingStrategy = ReportDFAndInfo()
- how and when to report results, seeReportingStrategy
name = "lomc!"
- name displayed in progress bar (viaProgressLogging
)metadata
- user-supplied metadata to be added to the reportdf
. Must be an iterable of pairs or aNamedTuple
, e.g.metadata = ("key1" => "value1", "key2" => "value2")
. All metadata is converted to strings.
Some metadata is automatically added to the report df
including Rimu.PACKAGE_VERSION
and data from state
.
Return values
lomc!
returns a named tuple with the following fields:
df
: aDataFrame
with all statistics being reported.state
: aQMCState
that can be used for continuations.
Example
julia> add = BoseFS(1,2,3);
julia> hamiltonian = HubbardReal1D(add);
julia> df1, state = lomc!(hamiltonian; targetwalkers=500, laststep=100);
julia> df2, _ = lomc!(state, df1; laststep=200, metadata=(;info="cont")); # Continuation run
julia> size(df1)
(100, 11)
julia> size(df2)
(200, 11)
julia> using DataFrames; metadata(df2, "info") # retrieve custom metadata
"cont"
julia> metadata(df2, "hamiltonian") # some metadata is automatically added
"HubbardReal1D(BoseFS{6,3}(1, 2, 3); u=1.0, t=1.0)"
Further keyword arguments and defaults:
τ_strat::TimeStepStrategy = ConstantTimeStep()
- adjust time step or not, seeTimeStepStrategy
s_strat::ShiftStrategy = DoubleLogUpdate(; targetwalkers, ζ = 0.08, ξ = ζ^2/4)
- how to update theshift
, seeShiftStrategy
.maxlength = 2 * s_strat.targetwalkers + 100
- upper limit on the length ofv
; when reached,lomc!
will abortparams::FciqmcRunStrategy = RunTillLastStep(laststep = 100, dτ = 0.01, shift = diagonal_element(ham, address)
- basic parameters of simulation state, seeFciqmcRunStrategy
. Parameter values are overridden by explicit keyword argumentslaststep
,dτ
,shift
; is mutated.wm
- working memory for re-use in subsequent calculations; is mutated.df = DataFrame()
- when called withAbstractHamiltonian
argument, aDataFrame
can be passed for merging with the reportdf
.
The default choice for the starting vector is v = default_starting_vector(; address, style, threading, initiator)
. See default_starting_vector
, PDVec
, DVec
, StochasticStyle
, and InitiatorRule
.
Rimu.post_step
— Functionpost_step(::PostStepStrategy, ::ReplicaState) -> kvpairs
Compute statistics after FCIQMC step. Should return a tuple of :key => value
pairs. This function is only called every reporting_interval
steps, as defined by the ReportingStrategy
.
See also PostStepStrategy
, ReportingStrategy
.
Rimu.refine_r_strat
— Methodrefine_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.
Rimu.replica_stats
— Functionreplica_stats(RS::ReplicaStrategy{N}, replicas::NTuple{N,ReplicaState}) -> (names, values)
Return the names and values of statistics related to N
replicas consistent with the ReplicaStrategy
RS
. names
should be a tuple of Symbol
s or String
s and values
should be a tuple of the same length. This function will be called every reporting_interval
steps from lomc!
, or once per time step if reporting_interval
is not defined.
Part of the ReplicaStrategy
interface. See also ReplicaState
.
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 at the very end of a step, after reporting_interval
steps. For example, it can be used to print some information to stdout
.
Rimu.report_metadata!
— Methodreport_metadata!(report::Report, key, value)
report_metadata!(report::Report, kvpairs)
Set metadata key
to value
in report
. key
and value
are converted to String
s. Alternatively, an iterable of key-value pairs or a NamedTuple
can be passed.
Throws an error if key
already exists.
See also get_metadata
, report!
, Report
.
Rimu.reporting_interval
— Methodreporting_interval(::ReportingStrategy)
Get the interval between steps for which non-essential statistics are reported. Defaults to 1 if chosen ReportingStrategy
does not specify an interval.
Rimu.single_particle_density
— Methodsingle_particle_density(dvec; component)
single_particle_density(add; component)
Compute the diagonal single particle density of vector dvec
or address add
. If the component
argument is given, only that component of the addresses is taken into account. The result is always normalized so that sum(result) ≈ num_particles(address)
.
Examples
julia> v = DVec(fs"|⋅↑⇅↓⋅⟩" => 1.0, fs"|↓↓⋅↑↑⟩" => 0.5)
DVec{FermiFS2C{2, 2, 5, 4, FermiFS{2, 5, BitString{5, 1, UInt8}}, FermiFS{2, 5, BitString{5, 1, UInt8}}},Float64} with 2 entries, style = IsDeterministic{Float64}()
fs"|↓↓⋅↑↑⟩" => 0.5
fs"|⋅↑⇅↓⋅⟩" => 1.0
julia> single_particle_density(v)
(0.2, 1.0, 1.6, 1.0, 0.2)
julia> single_particle_density(v; component=1)
(0.0, 1.6, 1.6, 0.4, 0.4)
See also
Rimu.smart_logger
— Methodsmart_logger(args...)
Enable terminal progress bar during interactive use (i.e. unless running on CI or HPC). Arguments are passed on to the logger. This is run once during Rimu
startup. Undo with default_logger
or by setting Base.global_logger()
.
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
Interfaces
StochasticStyles
Hamiltonians
BitStringAddresses
DictVectors
StatsTools
RMPI
See Module RMPI
Index
Rimu.DictVectors
Rimu.Hamiltonians
Rimu.Interfaces
Rimu.RMPI
Rimu.Rimu
Rimu.RimuIO
Rimu.StatsTools
Rimu.StochasticStyles
Rimu.PACKAGE_VERSION
Rimu.RMPI.mpi_root
Base.Matrix
Core.NamedTuple
MonteCarloMeasurements.Particles
Rimu.AllOverlaps
Rimu.BitStringAddresses.AbstractFockAddress
Rimu.BitStringAddresses.BitString
Rimu.BitStringAddresses.BoseFS
Rimu.BitStringAddresses.BoseFS2C
Rimu.BitStringAddresses.BoseFSIndex
Rimu.BitStringAddresses.CompositeFS
Rimu.BitStringAddresses.FermiFS
Rimu.BitStringAddresses.FermiFS2C
Rimu.BitStringAddresses.FermiFSIndex
Rimu.BitStringAddresses.OccupationNumberFS
Rimu.BitStringAddresses.OccupiedModeMap
Rimu.BitStringAddresses.OccupiedPairsMap
Rimu.BitStringAddresses.SingleComponentFockAddress
Rimu.BitStringAddresses.SortedParticleList
Rimu.ConstantTimeStep
Rimu.DictVectors.AbstractInitiatorValue
Rimu.DictVectors.AbstractProjector
Rimu.DictVectors.CoherentInitiator
Rimu.DictVectors.Communicator
Rimu.DictVectors.DVec
Rimu.DictVectors.FrozenDVec
Rimu.DictVectors.FrozenPDVec
Rimu.DictVectors.Initiator
Rimu.DictVectors.InitiatorDVec
Rimu.DictVectors.InitiatorRule
Rimu.DictVectors.InitiatorValue
Rimu.DictVectors.LocalPart
Rimu.DictVectors.MainSegmentIterator
Rimu.DictVectors.NonInitiator
Rimu.DictVectors.NonInitiatorValue
Rimu.DictVectors.Norm1ProjectorPPop
Rimu.DictVectors.Norm2Projector
Rimu.DictVectors.NormProjector
Rimu.DictVectors.NotDistributed
Rimu.DictVectors.PDVec
Rimu.DictVectors.PDWorkingMemory
Rimu.DictVectors.PDWorkingMemoryColumn
Rimu.DictVectors.PointToPoint
Rimu.DictVectors.SegmentedBuffer
Rimu.DictVectors.SimpleInitiator
Rimu.DictVectors.UniformProjector
Rimu.DontUpdate
Rimu.DoubleLogProjected
Rimu.DoubleLogSumUpdate
Rimu.DoubleLogUpdate
Rimu.DoubleLogUpdateAfterTargetWalkers
Rimu.FciqmcRunStrategy
Rimu.FirstOrderTransitionOperator
Rimu.Hamiltonians.AbstractOffdiagonals
Rimu.Hamiltonians.AxialAngularMomentumHO
Rimu.Hamiltonians.BasisSetRep
Rimu.Hamiltonians.BoseHubbardMom1D2C
Rimu.Hamiltonians.BoseHubbardReal1D2C
Rimu.Hamiltonians.DensityMatrixDiagonal
Rimu.Hamiltonians.ExtendedHubbardReal1D
Rimu.Hamiltonians.G2MomCorrelator
Rimu.Hamiltonians.G2RealCorrelator
Rimu.Hamiltonians.GuidingVectorSampling
Rimu.Hamiltonians.GutzwillerSampling
Rimu.Hamiltonians.HOCartesianCentralImpurity
Rimu.Hamiltonians.HOCartesianContactInteractions
Rimu.Hamiltonians.HOCartesianEnergyConservedPerDim
Rimu.Hamiltonians.HardwallBoundaries
Rimu.Hamiltonians.HubbardMom1D
Rimu.Hamiltonians.HubbardMom1DEP
Rimu.Hamiltonians.HubbardReal1D
Rimu.Hamiltonians.HubbardReal1DEP
Rimu.Hamiltonians.HubbardRealSpace
Rimu.Hamiltonians.LadderBoundaries
Rimu.Hamiltonians.LatticeGeometry
Rimu.Hamiltonians.MatrixHamiltonian
Rimu.Hamiltonians.Momentum
Rimu.Hamiltonians.Offdiagonals
Rimu.Hamiltonians.ParitySymmetry
Rimu.Hamiltonians.PeriodicBoundaries
Rimu.Hamiltonians.Stoquastic
Rimu.Hamiltonians.StringCorrelator
Rimu.Hamiltonians.SuperfluidCorrelator
Rimu.Hamiltonians.TimeReversalSymmetry
Rimu.Hamiltonians.Transcorrelated1D
Rimu.Interfaces.AbstractDVec
Rimu.Interfaces.AbstractHamiltonian
Rimu.Interfaces.CompressionStrategy
Rimu.Interfaces.LOStructure
Rimu.Interfaces.NoCompression
Rimu.Interfaces.StochasticStyle
Rimu.Interfaces.StyleUnknown
Rimu.LogUpdate
Rimu.LogUpdateAfterTargetWalkers
Rimu.MultiScalar
Rimu.NoStats
Rimu.PostStepStrategy
Rimu.ProjectedEnergy
Rimu.Projector
Rimu.QMCState
Rimu.RMPI.MPIAllToAll
Rimu.RMPI.MPIData
Rimu.RMPI.MPINoWalkerExchange
Rimu.RMPI.MPIOneSided
Rimu.RMPI.MPIPointToPoint
Rimu.ReplicaState
Rimu.ReplicaStrategy
Rimu.Report
Rimu.ReportDFAndInfo
Rimu.ReportToFile
Rimu.ReportingStrategy
Rimu.RunTillLastStep
Rimu.ShiftStrategy
Rimu.SignCoherence
Rimu.SingleParticleDensity
Rimu.StatsTools.BlockingResult
Rimu.StatsTools.RatioBlockingResult
Rimu.StochasticStyles.Bernoulli
Rimu.StochasticStyles.DynamicSemistochastic
Rimu.StochasticStyles.Exact
Rimu.StochasticStyles.IsDeterministic
Rimu.StochasticStyles.IsDynamicSemistochastic
Rimu.StochasticStyles.IsStochastic2Pop
Rimu.StochasticStyles.IsStochasticInteger
Rimu.StochasticStyles.IsStochasticWithThreshold
Rimu.StochasticStyles.SingleSpawn
Rimu.StochasticStyles.SpawningStrategy
Rimu.StochasticStyles.ThresholdCompression
Rimu.StochasticStyles.WithReplacement
Rimu.StochasticStyles.WithoutReplacement
Rimu.TimeStepStrategy
Rimu.Timer
Rimu.TripleLogUpdate
Rimu.WalkerLoneliness
Base.adjoint
LinearAlgebra.dot
Measurements.measurement
OrderedCollections.freeze
Rimu.BitStringAddresses.bose_hubbard_interaction
Rimu.BitStringAddresses.excitation
Rimu.BitStringAddresses.excitation
Rimu.BitStringAddresses.find_mode
Rimu.BitStringAddresses.find_occupied_mode
Rimu.BitStringAddresses.hopnextneighbour
Rimu.BitStringAddresses.near_uniform
Rimu.BitStringAddresses.num_components
Rimu.BitStringAddresses.num_modes
Rimu.BitStringAddresses.num_occupied_modes
Rimu.BitStringAddresses.num_particles
Rimu.BitStringAddresses.occupied_modes
Rimu.BitStringAddresses.onr
Rimu.BitStringAddresses.time_reverse
Rimu.DictVectors.collect_local!
Rimu.DictVectors.copy_to_local!
Rimu.DictVectors.dot_from_right
Rimu.DictVectors.from_initiator_value
Rimu.DictVectors.initiator_valtype
Rimu.DictVectors.is_distributed
Rimu.DictVectors.local_segments
Rimu.DictVectors.main_column
Rimu.DictVectors.merge_remote_reductions
Rimu.DictVectors.move_and_compress!
Rimu.DictVectors.mpi_comm
Rimu.DictVectors.mpi_comm
Rimu.DictVectors.mpi_rank
Rimu.DictVectors.mpi_rank
Rimu.DictVectors.mpi_recv_any!
Rimu.DictVectors.mpi_send
Rimu.DictVectors.mpi_size
Rimu.DictVectors.mpi_size
Rimu.DictVectors.num_columns
Rimu.DictVectors.num_rows
Rimu.DictVectors.perform_spawns!
Rimu.DictVectors.remote_segments
Rimu.DictVectors.replace_collections!
Rimu.DictVectors.synchronize_remote!
Rimu.DictVectors.synchronize_remote!
Rimu.DictVectors.target_segment
Rimu.DictVectors.to_initiator_value
Rimu.DictVectors.total_num_segments
Rimu.DictVectors.walkernumber
Rimu.Hamiltonians.build_basis
Rimu.Hamiltonians.check_address_type
Rimu.Hamiltonians.continuum_dispersion
Rimu.Hamiltonians.dimension
Rimu.Hamiltonians.fock_to_cart
Rimu.Hamiltonians.four_oscillator_integral_general
Rimu.Hamiltonians.get_all_blocks
Rimu.Hamiltonians.ho_delta_potential
Rimu.Hamiltonians.hubbard_dispersion
Rimu.Hamiltonians.log_abs_oscillator_zero
Rimu.Hamiltonians.momentum
Rimu.Hamiltonians.neighbour_site
Rimu.Hamiltonians.num_dimensions
Rimu.Hamiltonians.num_neighbours
Rimu.Hamiltonians.rayleigh_quotient
Rimu.Hamiltonians.shift_lattice
Rimu.Hamiltonians.shift_lattice_inv
Rimu.Interfaces.allowed_address_type
Rimu.Interfaces.apply_column!
Rimu.Interfaces.apply_operator!
Rimu.Interfaces.compress!
Rimu.Interfaces.default_style
Rimu.Interfaces.deposit!
Rimu.Interfaces.diagonal_element
Rimu.Interfaces.get_offdiagonal
Rimu.Interfaces.has_adjoint
Rimu.Interfaces.localpart
Rimu.Interfaces.num_offdiagonals
Rimu.Interfaces.offdiagonals
Rimu.Interfaces.random_offdiagonal
Rimu.Interfaces.sort_into_targets!
Rimu.Interfaces.starting_address
Rimu.Interfaces.step_stats
Rimu.Interfaces.storage
Rimu.Interfaces.working_memory
Rimu.RMPI.is_mpi_root
Rimu.RMPI.mpi_all_to_all
Rimu.RMPI.mpi_allprintln
Rimu.RMPI.mpi_barrier
Rimu.RMPI.mpi_combine_walkers!
Rimu.RMPI.mpi_no_exchange
Rimu.RMPI.mpi_one_sided
Rimu.RMPI.mpi_point_to_point
Rimu.RMPI.mpi_seed!
Rimu.RMPI.next_mpiID
Rimu.RMPI.targetrank
Rimu.RimuIO.load_df
Rimu.RimuIO.load_dvec
Rimu.RimuIO.save_df
Rimu.RimuIO.save_dvec
Rimu.StatsTools.autocovariance
Rimu.StatsTools.blocker
Rimu.StatsTools.blocking_analysis
Rimu.StatsTools.blocking_analysis_data
Rimu.StatsTools.blocks_with_m
Rimu.StatsTools.errs
Rimu.StatsTools.growth_estimator
Rimu.StatsTools.growth_estimator_analysis
Rimu.StatsTools.growth_witness
Rimu.StatsTools.growth_witness
Rimu.StatsTools.mean_and_se
Rimu.StatsTools.mixed_estimator
Rimu.StatsTools.mixed_estimator_analysis
Rimu.StatsTools.mtest
Rimu.StatsTools.particles
Rimu.StatsTools.particles
Rimu.StatsTools.projected_energy
Rimu.StatsTools.pseudo_cov
Rimu.StatsTools.ratio_estimators
Rimu.StatsTools.ratio_of_means
Rimu.StatsTools.rayleigh_replica_estimator
Rimu.StatsTools.rayleigh_replica_estimator_analysis
Rimu.StatsTools.replica_fidelity
Rimu.StatsTools.shift_estimator
Rimu.StatsTools.smoothen
Rimu.StatsTools.to_measurement
Rimu.StatsTools.val
Rimu.StatsTools.val_and_errs
Rimu.StatsTools.variational_energy_estimator
Rimu.StatsTools.w_exp
Rimu.StatsTools.w_lin
Rimu.StatsTools.x_by_y_linear
Rimu.StochasticStyles.diagonal_step!
Rimu.StochasticStyles.projected_deposit!
Rimu.StochasticStyles.spawn!
Rimu._n_walkers
Rimu.advance!
Rimu.all_overlaps
Rimu.check_transform
Rimu.default_logger
Rimu.default_starting_vector
Rimu.finalize_report!
Rimu.get_metadata
Rimu.lomc!
Rimu.post_step
Rimu.refine_r_strat
Rimu.replica_stats
Rimu.report!
Rimu.report!
Rimu.report_after_step
Rimu.report_metadata!
Rimu.reporting_interval
Rimu.single_particle_density
Rimu.smart_logger
Rimu.update_dτ
Rimu.update_shift
SparseArrays.sparse
Statistics.cov
Rimu.BitStringAddresses.@fs_str
Rimu.RMPI.@mpi_root