API
Rimu
Rimu.Rimu — ModuleRimuRandom integrators for many-body quantum systems
Welcome to Rimu version 0.15.0. Read the documentation online.
Rimu.PACKAGE_VERSION — ConstantRimu.PACKAGE_VERSIONConstant that contains the current VersionNumber of Rimu.
DataFrames.DataFrame — MethodDataFrame(report::Report)Convert the Report to a DataFrame. Metadata is added as metadata to the DataFrame.
Rimu.AllOverlaps — TypeAllOverlaps(n_replicas=2; operator=nothing, vecnorm=true, mixed_spectral_overlaps=false)
    <: ReplicaStrategy{n_replicas}Run n_replicas replicas and report overlaps between all pairs of replica vectors. If operator is not nothing, the overlap dot(r1, operator, r2) 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 r{i}s{k}_dot_r{j}s{k} for vector-vector overlaps, and r{i}s{k}_Op{m}_r{j}s{k} for operator overlaps, where i and j label the replicas, k labels the spectral state, and m labels the operators.
The r{i}s{k}_dot_r{j}s{k} overlap can be omitted with the flag vecnorm=false.
By default, overlaps of different spectral states are omitted. To include overlaps of different spectral states r{i}s{k}_dot_r{j}s{l} and r{i}s{k}_Op{m}_r{j}s{l}, use the flag mixed_spectral_overlaps=true.
See ProjectorMonteCarloProblem, ReplicaStrategy and AbstractOperator (for an interface for implementing operators).
Transformed Hamiltonians
If a transformed Hamiltonian G has been passed to ProjectorMonteCarloProblem, an inverse transformation is applied to the operators in AllOverlaps. Additionally, an operator representing the inverse transform applied to the identity operator is added to the list of operators. Passing transform to AllOverlaps is deprecated.
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.)
\[ Ĝ = f Ĥ f⁻¹.\]
The expectation value of an operator $Â$ is
\[ ⟨Â⟩ = ⟨ψ| Â |ψ⟩ = \frac{⟨ϕ| f⁻¹ Â f⁻¹ |ϕ⟩}{⟨ϕ| f⁻² |ϕ⟩}\]
where
\[ |ϕ⟩ = f |ψ⟩\]
is the (right) eigenvector of $Ĝ$ and $|ψ⟩$ is an eigenvector of $Ĥ$.
For an $m$-tuple of input operators $(\hat{A}_1, ..., \hat{A}_m)$, overlaps of $⟨ϕ| f⁻¹ Â f⁻¹ |ϕ⟩$ are reported as r{i}s{k}_Op{m}_r{j}s{k}. The correct vector-vector overlap $⟨ϕ| f⁻² |ϕ⟩$ is reported last as r{i}s{k}_Op{m+1}_r{j}s{k}. This is in addition to the bare vector-vector overlap $⟨ϕ|ϕ⟩$ that is reported as r{i}s{k}_dot_r{j}s{k}.
Rimu.ConstantTimeStep — TypeConstantTimeStep <: TimeStepStrategyKeep time_step constant.
Rimu.DefaultShiftParameters — TypeDefaultShiftParametersDefault mutable struct for storing the shift parameters.
Rimu.DontUpdate — TypeDontUpdate(; target_walkers = 1_000) <: ShiftStrategyDon't update the shift.  Return when target_walkers 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)\]
Note that adjusting the keyword maxlength in ProjectorMonteCarloProblem is advised as the default may not be appropriate.
Rimu.DoubleLogSumUpdate — TypeDoubleLogSumUpdate(; target_walkers = 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(; target_walkers = 1_000, ζ = 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(target_walkers = 1_000, ζ = 0.08, ξ = ζ^2/4) <: ShiftStrategyStrategy for updating the shift: After target_walkers is reached, update the shift according to the log formula with damping parameter ζ and ξ.
See DoubleLogUpdate, ShiftStrategy, ProjectorMonteCarloProblem.
Rimu.FCIQMC — TypeFCIQMC(; kwargs...) <: PMCAlgorithmAlgorithm for the full configuration interaction quantum Monte Carlo (FCIQMC) method. The default algorithm for ProjectorMonteCarloProblem.
Keyword arguments and defaults:
- shift_strategy = DoubleLogUpdate(; targetwalkers = 1_000, ζ = 0.08, ξ = ζ^2/4): How to update the- shift.
- time_step_strategy = ConstantTimeStep(): Adjust time step or not.
See also ProjectorMonteCarloProblem, ShiftStrategy, TimeStepStrategy, DoubleLogUpdate, ConstantTimeStep.
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:
The use of this strategy is deprecated. Pass the relevant arguments directly to ProjectorMonteCarloProblem or to lomc!() instead.
Rimu.FirstOrderTransitionOperator — TypeFirstOrderTransitionOperator(hamiltonian, shift, time_step) <: AbstractHamiltonian
FirstOrderTransitionOperator(sp::DefaultShiftParameters, hamiltonian)First order transition operator
\[𝐓 = 1 + dτ(S - 𝐇)\]
where $𝐇$ is the hamiltonian, $dτ$ the time_step 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.GramSchmidt — TypeGramSchmidt(S; orthogonalization_interval = 1) <: SpectralStrategy{S}Use the Gram-Schmidt procedure to orthogonalize the excited states. A total of S spectral states are used in the simulation, and they are orthogonalized every  orthogonalization_interval steps.
Use with the keyword argument spectral_strategy in ProjectorMonteCarloProblem.
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(target_walkers = 1_000, ζ = 0.08) <: ShiftStrategyStrategy for updating the shift: After target_walkers is reached, update the shift according to the log formula with damping parameter ζ.
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)Note that only MultiScalars 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 ProjectorMonteCarloProblem.
Rimu.PMCAlgorithm — TypePMCAlgorithmAbstract type for projector Monte Carlo algorithms.
Rimu.PMCSimulation — TypePMCSimulationHolds the state and the results of a projector quantum Monte Carlo (PMC) simulation. Is returned by init(::ProjectorMonteCarloProblem) and solved with solve!(::PMCSimulation).
Obtain the results of a simulation sm as a DataFrame with DataFrame(sm).
Fields
- problem::ProjectorMonteCarloProblem: The problem that was solved
- state::Rimu.ReplicaState: The current state of the simulation
- report::Rimu.Report: The report of the simulation
- modified::Bool: Whether the simulation has been modified
- aborted::Bool: Whether the simulation has been aborted
- success::Bool: Whether the simulation has been completed successfully
- message::String: A message about the simulation status
- elapsed_time::Float64: The time elapsed during the simulation
See also state_vectors, ProjectorMonteCarloProblem, init, solve!.
Rimu.PostStepStrategy — TypePostStepStrategySubtypes of PostStepStrategy can be used to perform arbitrary computation on a single state after an FCIQMC step is finished and report the results.
Implemented strategies:
Note: A tuple of multiple strategies can be passed to ProjectorMonteCarloProblem. In that case, all reported column names must be distinct.
Interface:
A subtype of this type must implement post_step_action(::PostStepStrategy, ::SingleState, step::Int).
Rimu.ProjectedEnergy — TypeProjectedEnergy(hamiltonian, projector; hproj=:hproj, 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 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) <: PostStepStrategyAfter each step, compute dot(projector, dv) and report it in the DataFrame under name. projector can be an AbstractDVec, or an AbstractProjector.
Rimu.ProjectorMonteCarloProblem — TypeProjectorMonteCarloProblem(hamiltonian::AbstractHamiltonian; kwargs...)Defines a problem to be solved by projector quantum Monte Carlo (QMC) methods, such as the the FCIQMC algorithm.
Common keyword arguments and defaults:
- time_step = 0.01: Initial time step size.
- last_step = 100: Controls the number of steps.
- target_walkers = 1_000: Target for the 1-norm of the coefficient vector.
- start_at = starting_address(hamiltonian): Define the initial state vector(s). An $r × s$ matrix of state vectors can be passed where $r$ is the number of replicas and $s$ the number of spectral states. See also- default_starting_vector.
- style = IsDynamicSemistochastic(): The- StochasticStyleof the simulation.
- initiator = false: Whether to use initiators. Can be- true,- false, or a valid- InitiatorRule.
- threading: Default is to use multithreading and/or MPI if available. Set to- trueto force- PDVecfor the starting vector,- falsefor serial computation; may be overridden by- start_at.
- reporting_strategy = ReportDFAndInfo(): How and when to report results, see- ReportingStrategy.
- post_step_strategy = (): Extract observables (e.g.- ProjectedEnergy), see- PostStepStrategy.
- n_replicas = 1: Number of synchronised independent simulations.
- replica_strategy = NoStats(n_replicas): Which results to report from replica simulations, see- ReplicaStrategy.
- n_spectral = 1: Number of targeted spectral states. Set- n_spectral > 1to find excited states.
- spectral_strategy = GramSchmidt(n_spectral): The- SpectralStrategyused for orthogonalizing spectral states.
Example
julia> hamiltonian = HubbardReal1D(BoseFS(1,2,3));
julia> problem = ProjectorMonteCarloProblem(hamiltonian; target_walkers = 500, last_step = 100);
julia> simulation = solve(problem);
julia> simulation.success[]
true
julia> size(DataFrame(simulation))
(100, 9)Further keyword arguments:
- starting_step = 1: Starting step of the simulation.
- wall_time = Inf: Maximum time allowed for the simulation.
- simulation_plan = SimulationPlan(; starting_step, last_step, wall_time): Defines the duration of the simulation. Takes precedence over- last_stepand- wall_time.
- ζ = 0.08: Damping parameter for the shift update.
- ξ = ζ^2/4: Forcing parameter for the shift update.
- shift_strategy = DoubleLogUpdate(; target_walkers, ζ, ξ): How to update the- shift, see- ShiftStrategy.
- time_step_strategy = ConstantTimeStep(): Adjust time step or not, see- TimeStepStrategy.
- algorithm = FCIQMC(; shift_strategy, time_step_strategy): The algorithm to use. Currenlty only- FCIQMCis implemented.
- shift: Initial shift value or collection of shift values. Determined by default from the Hamiltonian and the starting vectors.
- initial_shift_parameters: Initial shift parameters or collection of initial shift parameters. Overrides- shiftif provided.
- max_length = 2 * target_walkers + 100: Maximum length of the vectors.
- display_name = "PMCSimulation": Name displayed in progress bar (via- ProgressLogging).
- metadata: User-supplied metadata to be added to the report. Must be an iterable of pairs or a- NamedTuple, e.g.- metadata = ("key1" => "value1", "key2" => "value2"). All metadata is converted to strings.
- random_seed = true: Provide and store a seed for the random number generator. If set to- true, a new random seed is generated from- RandomDevice(). If set to number, this number is used as the seed. This seed is used by- solve(and- init) to re-seed the default random number generator (consistently on each MPI rank) such that- solveing the same- ProjectorMonteCarloProblemtwice will yield identical results. If set to- false, no seed is used and consecutive random numbers are used.
- minimum_size = 2*num_spectral_states(spectral_strategy): The minimum size of the basis used to construct starting vectors for simulations of spectral states, if- start_atis not provided.
Rimu.ReplicaState — TypeReplicaState <: AbstractMatrix{SingleState}Holds information about multiple replicas of SpectralStates. Indexing the ReplicaState state[i, j] returns a SingleState from the ith replica and jth spectral state.
Fields
- spectral_states: Tuple of- SpectralStates
- max_length::Ref{Int}: Maximum length of the simulation
- step::Ref{Int}: Current step of the simulation
- simulation_plan: Simulation plan
- reporting_strategy: Reporting strategy
- post_step_strategy: Post-step strategy
- replica_strategy: Replica strategy
See also ReplicaStrategy, Rimu.SpectralState, Rimu.SingleState, Rimu.PMCSimulation.
Rimu.ReplicaStrategy — TypeReplicaStrategy{N}Supertype for strategies that can be passed to ProjectorMonteCarloProblem 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 of- Strings or- Symbolsof names for replica statistics and a tuple of the values. These will be reported to the- DataFramereturned by- ProjectorMonteCarloProblem.
Rimu.Report — TypeReport()
Report(df::DataFrame)Internal structure that holds the temporary reported values as well as metadata. It can be converted to a DataFrame with DataFrame(report::Report).
See report!, report_metadata!, get_metadata.
Rimu.ReportDFAndInfo — TypeReportDFAndInfo(; reporting_interval=1, info_interval=100, io=stdout, writeinfo=false) <: ReportingStrategyThe default ReportingStrategy for ProjectorMonteCarloProblem. Report every reporting_intervalth step to a DataFrame and write info message to io every info_intervalth 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().
See also ProjectorMonteCarloProblem, ReportToFile.
Rimu.ReportToFile — TypeReportToFile(; kwargs...) <: ReportingStrategyReportingStrategy for ProjectorMonteCarloProblem 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. A- DataFrameof this size is collected in memory and written to disk. When saving, an info message is also printed to- io.
- 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 empty- DataFrameis returned.
- io = stdout: The- IOto print messages to. Set to- devnullif you don't want to see messages printed out.
- compress = :zstd: compression algorithm to use. Can be- :zstd,- :lz4or- nothing.
See also load_df, save_df, ReportDFAndInfo, and ProjectorMonteCarloProblem.
Rimu.ReportingStrategy — TypeReportingStrategyAbstract type for strategies for reporting data during a simulation of a ProjectorMonteCarloProblem.
Implemented strategies:
Extended help
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 lomc!() for a fixed number of time steps. For alternative strategies, see FciqmcRunStrategy.
The use of this strategy is deprecated. Pass the relevant arguments directly to ProjectorMonteCarloProblem or to lomc!() instead.
Rimu.ShiftStrategy — TypeShiftStrategyAbstract type for defining the strategy for controlling the norm, potentially by updating the shift. Passed as a parameter to ProjectorMonteCarloProblem or to FCIQMC.
Implemented strategies:
- DontUpdate
- DoubleLogUpdate- default in- ProjectorMonteCarloProblem()
- LogUpdate
- LogUpdateAfterTargetWalkers- FCIQMC standard
- DoubleLogUpdateAfterTargetWalkers
Extended help
Internally To implement a custom strategy, define a new subtype of Rimu.ShiftStrategy and implement methods for:
- Rimu.update_shift_parameters!- to update the- shift_parameters
- Rimu.initialise_shift_parameters- (optional) to initialise and construct a custom implementation of the- shift_parameters. The default implementation is- Rimu.DefaultShiftParameters.
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.SimulationPlan — TypeSimulationPlan(; starting_step = 1, last_step = 100, wall_time = Inf)Defines the duration of the simulation. The simulation ends when the last_step is reached or the wall_time is exceeded.
Rimu.SingleParticleDensity — TypeSingleParticleDensity(; save_every=1, component) <: PostStepStrategyPostStepStrategy  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.SingleState — TypeSingleState(hamiltonian, algorithm, v, wm, pnorm, params, id)Struct that holds a single state vector and all information needed for an independent run of the algorithm. Can be advanced a step forward with Rimu.advance!.
Fields
- hamiltonian: Hamiltonian
- algorithm: Algorithm
- v: Vector
- pv: Previous vector
- wm: Working memory
- shift_parameters: Shift parameters
- id::String: id is appended to column names
See also SpectralStrategy, ReplicaStrategy, Rimu.SpectralState, Rimu.ReplicaState, Rimu.replica_stats, Rimu.PMCSimulation.
Rimu.SpectralState — TypeSpectralState <: AbstractVector{SingleState}Holds one or several Rimu.SingleStates representing the ground state and excited states of a single replica. Indexing the SpectralState state[i] returns the ith SingleState.
Fields
- single_states: Tuple of- SingleStates
- spectral_strategy: Strategy for computing the spectral states
- id::String: Identifies the replica
See also SpectralStrategy, Rimu.ReplicaState, Rimu.SingleState, Rimu.PMCSimulation.
Rimu.SpectralStrategy — TypeSpectralStrategy{S}Abstract type for spectral strategies. The spectral strategy is used to control the number of spectral states used in the simulation.
Implemented Strategies
- GramSchmidt: Orthogonalize the spectral states using the Gram-Schmidt procedure.
Rimu.StateVectors — TypeStateVectors <: AbstractMatrix{V}Represents a matrix of configuration vectors from the state. Construct this object with state_vectors.
Rimu.TimeStepStrategy — TypeTimeStepStrategyAbstract type for strategies for updating the time step with update_time_step(). Implemented strategies:
Rimu.Timer — TypeTimer <: PostStepStrategyRecord current time after every step. See Base.Libc.time for information on what time is recorded.
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.
CommonSolve.init — Methodinit(problem::ProjectorMonteCarloProblem; copy_vectors=true)::PMCSimulationInitialise a Rimu.PMCSimulation.
See also ProjectorMonteCarloProblem, solve!, solve, step!, Rimu.PMCSimulation.
CommonSolve.solve — Functionsolve(::ProjectorMonteCarloProblem)::PMCSimulationInitialize and solve a ProjectorMonteCarloProblem until the last step is completed or the wall time limit is reached.
See also init, solve!, step!, Rimu.PMCSimulation, and solve(::ExactDiagonalizationProblem).
CommonSolve.solve! — Methodsolve!(sm::PMCSimulation; kwargs...)::PMCSimulationSolve a Rimu.PMCSimulation until the last step is completed or the wall time limit is reached.
To continue a previously completed simulation, set a new last_step or wall_time using the keyword arguments. Optionally, changes can be made to the replica_strategy, the post_step_strategy, or the reporting_strategy.
Optional keyword arguments:
- last_step = nothing: Set the last step to a new value and continue the simulation.
- wall_time = nothing: Set the allowed wall time to a new value and continue the simulation.
- reset_time = false: Reset the- elapsed_timecounter and continue the simulation.
- empty_report = false: Empty the report before continuing the simulation.
- replica_strategy = nothing: Change the replica strategy. Requires the number of replicas to match the number of replicas in the simulation- sm. Implies- empty_report = true.
- post_step_strategy = nothing: Change the post-step strategy. Implies- empty_report = true.
- reporting_strategy = nothing: Change the reporting strategy. Implies- empty_report = true.
- metadata = nothing: Add metadata to the report.
See also ProjectorMonteCarloProblem, init, solve, step!, Rimu.PMCSimulation.
CommonSolve.step! — Methodstep!(sm::PMCSimulation)::PMCSimulationAdvance the simulation by one step.
Calling solve! will advance the simulation until the last step or the wall time is exceeded. When completing the simulation without calling solve!, the simulation report needs to be finalised by calling Rimu.finalize_report!.
See also ProjectorMonteCarloProblem, init, solve!, solve, Rimu.PMCSimulation.
Rimu.Interfaces.num_overlaps — Methodnum_overlaps(state_or_DataFrame)Return the number of overlaps between replicas that are being reported. Only counts  overlaps between replicas of the same spectral state, even if AllOverlaps is used with  mixed_spectral_overlaps=true.
The return value is 0 if no overlaps are being reported, and N*(N-1)÷2 otherwise, where N is the value returned by num_replicas.
See also ProjectorMonteCarloProblem, AllOverlaps, num_spectral_states.
Rimu.Interfaces.num_replicas — Methodnum_replicas(state_or_strategy_or_DataFrame)Return the number of replicas used in the simulation. With multiple spectral states, only reports the number of replicas per spectral state.
See also ProjectorMonteCarloProblem, AllOverlaps, num_spectral_states.
Rimu.Interfaces.num_spectral_states — Methodnum_spectral_states(state_or_strategy_or_DataFrame)Return the number of spectral states used in the simulation.
See also ProjectorMonteCarloProblem, GramSchmidt, num_replicas.
Rimu.advance! — Methodadvance!(algorithm::PMCAlgorithm, report::Report, state::ReplicaState, s_state::SingleState)Advance the s_state by one step according to the algorithm. 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, mixed_spectral_overlaps=false)Get all overlaps between vectors and operators.  The flag vecnorm can disable the vector-vector overlap r{i}s{k}_dot_r{j}s{k}.
Rimu.clean_and_warn_if_others_present — Methodclean_and_warn_if_others_present(nt::NamedTuple{names}, keys) where {names}Remove keys from a NamedTuple that are not in keys and issue a warning if they are present.
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=IsDynamicSemistochastic(),
    initiator=NonInitiator(),
    threading=nothing,
    population=10
)Return a default starting vector for ProjectorMonteCarloProblem. The default choice for the starting vector is
v = PDVec(address => population; style, initiator)if threading is available, or otherwise
v = DVec(address => population; style)if initiator == NonInitiator(), and
v = InitiatorDVec(address => population; style, initiator)if not. See PDVec, DVec, InitiatorDVec, StochasticStyle, and InitiatorRule.
Rimu.delete_and_warn_if_present — Methoddelete_and_warn_if_present(nt::NamedTuple, keys)Delete keys from a NamedTuple and issue a warning if they are present. This is useful for removing unused keyword arguments.
Rimu.finalize_report! — Methodfinalize_report!(::ReportingStrategy, report)Finalize the report. This function is called after all steps in solve! 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.initialise_shift_parameters — Functioninitialise_shift_parameters(s::ShiftStrategy, shift, norm, time_step, counter=0, shift_mode=false)Initiatlise a struct to store the shift parameters.
See ShiftStrategy, update_shift_parameters!, DefaultShiftParameters.
Rimu.is_mpi_root — Functionis_mpi_root(root = mpi_root)Returns true if called from the root rank
Rimu.lomc! — Methodlomc!(ham::AbstractHamiltonian, [v]; kwargs...) -> df, state
lomc!(state::ReplicaState, [df]; kwargs...) -> df, stateLinear 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 ReplicaState 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 default- vand- shift.
- style = IsStochasticInteger()- set- StochasticStylefor default- v; unused if- vis specified.
- initiator = NonInitiator()- set- InitiatorRulefor default- v; unused if- vis specified.
- threading- default is to use multithreading and MPI if multiple threads are available. Set to- trueto force- PDVecfor the starting vector,- falsefor serial computation; unused if- vis specified.
- shift = diagonal_element(ham, address)- initial value of shift.
- post_step_strategy::NTuple{N,<:PostStepStrategy} = ()- extract observables (e.g.- ProjectedEnergy), see- PostStepStrategy. (Deprecated:- post_stepis accepted as an alias for- post_step_strategy.)
- replica_strategy::ReplicaStrategy = NoStats(1)- run several synchronised simulations, see- ReplicaStrategy. (Deprecated:- replicais accepted as an alias for- replica_strategy.)
- reporting_strategy::ReportingStrategy = ReportDFAndInfo()- how and when to report results, see- ReportingStrategy. (Deprecated:- r_stratis accepted as an alias for- reporting_strategy.)
- name = "lomc!"- name displayed in progress bar (via- ProgressLogging)
- metadata- user-supplied metadata to be added to the report- df. Must be an iterable of pairs or a- NamedTuple, 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: a- DataFramewith all statistics being reported.
- state: a- ReplicaStatethat can be used for continuations.
Example
julia> address = BoseFS(1,2,3);
julia> hamiltonian = HubbardReal1D(address);
julia> df1, state = @suppress lomc!(hamiltonian; targetwalkers=500, laststep=100);
julia> df2, _ = @suppress lomc!(state, df1; laststep=200, metadata=(;info="cont")); # Continuation run
julia> size(df1)
(100, 9)
julia> size(df2)
(200, 9)
julia> using DataFrames; metadata(df2, "info") # retrieve custom metadata
"cont"
julia> metadata(df2, "hamiltonian") # some metadata is automatically added
"HubbardReal1D(fs\"|1 2 3⟩\"; u=1.0, t=1.0)"Further keyword arguments and defaults:
- τ_strat::TimeStepStrategy = ConstantTimeStep()- adjust time step or not, see- TimeStepStrategy
- s_strat::ShiftStrategy = DoubleLogUpdate(; target_walkers=targetwalkers, ζ = 0.08, ξ = ζ^2/4)- how to update the- shift, see- ShiftStrategy.
- maxlength = 2 * s_strat.target_walkers + 100- upper limit on the length of- v; when reached,- lomc!will abort
- wm- working memory for re-use in subsequent calculations; is mutated.
- df = DataFrame()- when called with- AbstractHamiltonianargument, a- DataFramecan be passed for merging with the report- df.
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.
The use of this lomc! is deprecated. Use ProjectorMonteCarloProblem and solve instead.
Rimu.mpi_allprintln — Methodmpi_allprintln(args...)Print a message to stdout from each rank separately, in order. MPI synchronizing.
Rimu.mpi_barrier — Functionmpi_barrier(comm = mpi_comm())The MPI barrier with optional argument. MPI syncronizing.
Rimu.mpi_rank — Functionmpi_rank(comm = mpi_comm())Return the current MPI rank.
Rimu.mpi_seed! — Functionmpi_seed!(seed = rand(Random.RandomDevice(), UInt))Re-seed the random number generators in an MPI-safe way. If seed is provided, the random numbers from rand 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.mpi_size — Functionmpi_size(comm = mpi_comm())Size of MPI communicator.
Rimu.post_step_action — Functionpost_step_action(::PostStepStrategy, ::SingleState, step) -> kvpairsCompute 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_reporting_strategy — Methodrefine_reporting_strategy(reporting_strategy::ReportingStrategy) -> reporting_strategyInitialize 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.replace_keys — Methodreplace_keys(nt::NamedTuple, (:old1 => :new1, :old2 => :new2, ...))Replace keys in a NamedTuple with new keys. This is useful for renaming fields in a NamedTuple. Ignores keys that are not present in the NamedTuple.
Rimu.replica_stats — Functionreplica_stats(RS::ReplicaStrategy{N}, spectral_states::NTuple{N,SingleState}) -> (names, values)Return the names and values of statistics related to N replica states consistent with the ReplicaStrategy RS. names should be a tuple of Symbols or Strings and values should be a tuple of the same length. This function will be called every reporting_interval steps from ProjectorMonteCarloProblem, or once per time step if reporting_interval is not defined.
Part of the ReplicaStrategy interface. See also SingleState.
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 ProjectorMonteCarloProblem 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::Report, df::DataFrame)Convert the DataFrame df to a Report. This function does not copy the data.
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) -> reportThis function is called at the very end of a step, after Rimu.reporting_interval steps. It may modify the report.
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 Strings. Alternatively, an iterable of key-value pairs or a NamedTuple can be passed.
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.set_up_initial_shift_parameters — Methodset_up_initial_shift_parameters(
    algorithm::FCIQMC, hamiltonian, starting_vectors, shift, time_step, initial_shift_parameters
)Set up the initial shift parameters for the FCIQMC algorithm.
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"|⋅↑⇅↓⋅⟩" => 1.0
  fs"|↓↓⋅↑↑⟩" => 0.5
julia> single_particle_density(v)
(0.2, 1.0, 1.6, 1.0, 0.2)
julia> single_particle_density(v; component=1)
(0.0, 0.8, 0.8, 0.2, 0.2)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.split_keys — Methodsplit_keys(nt::NamedTuple, selected_keys...) -> extracted, restFind keys from selected_keys that have values in nt. Return them as extracted and return nt with those keys deleted as rest.
Example
julia> nt = (;a=1, b=2, d=3);
julia> Rimu.split_keys(nt, :a, :c)
((a = 1,), (b = 2, d = 3))Rimu.state_vectors — Methodstate_vectors(state::ReplicaState)
state_vectors(sim::PMCSimulation)Return an r×s AbstractMatrix of configuration vectors from the state, or the result of solve(::ProjectorMonteCarloProblem). The vectors can be accessed by indexing the resulting collection, where the row index corresponds to the replica index and the column index corresponds to the spectral state index.
See also ProjectorMonteCarloProblem, Rimu.PMCSimulation, Rimu.SingleState, Rimu.ReplicaState, Rimu.SpectralState.
Rimu.update_shift_parameters! — Functionupdate_shift_parameters!(
    s <: ShiftStrategy,
    shift_parameters,
    tnorm,
    v_new,
    v_old,
    step,
    report
) -> shift_stats, proceedUpdate the shift_parameters according to strategy s. See ShiftStrategy. Returns a named tuple of the shift statistics and a boolean proceed indicating whether the simulation should proceed.
Rimu.update_time_step — Methodupdate_time_step(s<:TimeStepStrategy, time_step, tnorm) -> new_time_stepUpdate the time step according to the strategy s.
Rimu.@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 onlyReexported Submodules
ExactDiagonalization
Interfaces
StochasticStyles
Hamiltonians
BitStringAddresses
DictVectors
StatsTools
Index
- Rimu.BitStringAddresses
- Rimu.DictVectors
- Rimu.ExactDiagonalization
- Rimu.Hamiltonians
- Rimu.InterfaceTests
- Rimu.Interfaces
- Rimu.Rimu
- Rimu.RimuIO
- Rimu.StatsTools
- Rimu.StochasticStyles
- Rimu.PACKAGE_VERSION
- Base.Matrix
- Core.NamedTuple
- DataFrames.DataFrame
- LinearMaps.LinearMap
- MonteCarloMeasurements.Particles
- Rimu.AllOverlaps
- Rimu.BitStringAddresses.BitString
- Rimu.BitStringAddresses.BoseFS
- Rimu.BitStringAddresses.BoseFSIndex
- Rimu.BitStringAddresses.CompositeFS
- Rimu.BitStringAddresses.FermiFS
- Rimu.BitStringAddresses.FermiFS2C
- Rimu.BitStringAddresses.FermiFS2CModes
- Rimu.BitStringAddresses.FermiFSIndex
- Rimu.BitStringAddresses.ModeMap
- Rimu.BitStringAddresses.OccupationNumberFS
- Rimu.BitStringAddresses.OccupiedPairsMap
- Rimu.BitStringAddresses.SingleComponentFockAddress
- Rimu.BitStringAddresses.SortedParticleList
- Rimu.ConstantTimeStep
- Rimu.DefaultShiftParameters
- Rimu.DictVectors.AbstractInitiatorValue
- Rimu.DictVectors.AbstractProjector
- Rimu.DictVectors.AllToAll
- Rimu.DictVectors.CoherentInitiator
- Rimu.DictVectors.Communicator
- Rimu.DictVectors.DVec
- Rimu.DictVectors.FirstColumnIterator
- Rimu.DictVectors.FrozenDVec
- Rimu.DictVectors.FrozenPDVec
- Rimu.DictVectors.Initiator
- Rimu.DictVectors.InitiatorDVec
- Rimu.DictVectors.InitiatorRule
- Rimu.DictVectors.InitiatorValue
- Rimu.DictVectors.LocalPart
- Rimu.DictVectors.NestedSegmentedBuffer
- 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.ExactDiagonalization.ArpackSolver
- Rimu.ExactDiagonalization.BasisSetRepresentation
- Rimu.ExactDiagonalization.ExactDiagonalizationProblem
- Rimu.ExactDiagonalization.KrylovKitSolver
- Rimu.ExactDiagonalization.LOBPCGSolver
- Rimu.ExactDiagonalization.LinearAlgebraSolver
- Rimu.ExactDiagonalization.OperatorAsMap
- Rimu.FCIQMC
- Rimu.FciqmcRunStrategy
- Rimu.FirstOrderTransitionOperator
- Rimu.GramSchmidt
- Rimu.Hamiltonians.AbstractOffdiagonals
- Rimu.Hamiltonians.AxialAngularMomentumHO
- Rimu.Hamiltonians.CubicGrid
- Rimu.Hamiltonians.DensityMatrixDiagonal
- Rimu.Hamiltonians.Directions
- Rimu.Hamiltonians.Displacements
- Rimu.Hamiltonians.ExtendedHubbardMom1D
- Rimu.Hamiltonians.ExtendedHubbardReal1D
- Rimu.Hamiltonians.FroehlichPolaron
- Rimu.Hamiltonians.G2MomCorrelator
- Rimu.Hamiltonians.G2RealCorrelator
- Rimu.Hamiltonians.G2RealSpace
- Rimu.Hamiltonians.GuidingVectorSampling
- Rimu.Hamiltonians.GutzwillerSampling
- Rimu.Hamiltonians.HOCartesianCentralImpurity
- Rimu.Hamiltonians.HOCartesianContactInteractions
- Rimu.Hamiltonians.HOCartesianEnergyConservedPerDim
- Rimu.Hamiltonians.HamiltonianProduct
- Rimu.Hamiltonians.HamiltonianSum
- Rimu.Hamiltonians.HubbardMom1D
- Rimu.Hamiltonians.HubbardMom1DEP
- Rimu.Hamiltonians.HubbardReal1D
- Rimu.Hamiltonians.HubbardReal1DEP
- Rimu.Hamiltonians.HubbardRealSpace
- Rimu.Hamiltonians.IdentityOperator
- Rimu.Hamiltonians.MatrixHamiltonian
- Rimu.Hamiltonians.ModifiedHamiltonian
- Rimu.Hamiltonians.MolecularHamiltonian
- Rimu.Hamiltonians.MolecularHamiltonianOffDiagonals
- Rimu.Hamiltonians.MolecularHamiltonianOffDiagonalsIteratorState
- Rimu.Hamiltonians.Momentum
- Rimu.Hamiltonians.Offdiagonals
- Rimu.Hamiltonians.ParitySymmetry
- Rimu.Hamiltonians.ParticleNumberOperator
- Rimu.Hamiltonians.ReducedDensityMatrix
- Rimu.Hamiltonians.ScaledHamiltonian
- Rimu.Hamiltonians.SingleParticleExcitation
- Rimu.Hamiltonians.Stoquastic
- Rimu.Hamiltonians.StringCorrelator
- Rimu.Hamiltonians.SuperfluidCorrelator
- Rimu.Hamiltonians.TimeReversalSymmetry
- Rimu.Hamiltonians.Transcorrelated1D
- Rimu.Hamiltonians.TransformUndoer
- Rimu.Hamiltonians.TwoParticleExcitation
- Rimu.Interfaces.AbstractDVec
- Rimu.Interfaces.AbstractFockAddress
- Rimu.Interfaces.AbstractHamiltonian
- Rimu.Interfaces.AbstractObservable
- Rimu.Interfaces.AbstractOperator
- Rimu.Interfaces.AbstractOperatorColumn
- Rimu.Interfaces.CompressionStrategy
- Rimu.Interfaces.LOStructure
- Rimu.Interfaces.NoCompression
- Rimu.Interfaces.OffdiagonalsOperatorColumn
- Rimu.Interfaces.StochasticStyle
- Rimu.Interfaces.StyleUnknown
- Rimu.LogUpdate
- Rimu.LogUpdateAfterTargetWalkers
- Rimu.MultiScalar
- Rimu.NoStats
- Rimu.PMCAlgorithm
- Rimu.PMCSimulation
- Rimu.PostStepStrategy
- Rimu.ProjectedEnergy
- Rimu.Projector
- Rimu.ProjectorMonteCarloProblem
- Rimu.ReplicaState
- Rimu.ReplicaStrategy
- Rimu.Report
- Rimu.ReportDFAndInfo
- Rimu.ReportToFile
- Rimu.ReportingStrategy
- Rimu.RimuIO.DVecAsTable
- Rimu.RimuIO.PDVecAsTable
- Rimu.RunTillLastStep
- Rimu.ShiftStrategy
- Rimu.SignCoherence
- Rimu.SimulationPlan
- Rimu.SingleParticleDensity
- Rimu.SingleState
- Rimu.SpectralState
- Rimu.SpectralStrategy
- Rimu.StateVectors
- 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.WalkerLoneliness
- Base.adjoint
- Base.eltype
- Base.mapreduce
- CommonSolve.init
- CommonSolve.init
- CommonSolve.solve
- CommonSolve.solve
- CommonSolve.solve!
- CommonSolve.step!
- LinearAlgebra.dot
- LinearAlgebra.mul!
- Measurements.measurement
- OrderedCollections.freeze
- Rimu.BitStringAddresses.bose_hubbard_interaction
- Rimu.BitStringAddresses.each_mode
- Rimu.BitStringAddresses.excitation
- Rimu.BitStringAddresses.excitation
- Rimu.BitStringAddresses.find_mode
- Rimu.BitStringAddresses.find_occupied_mode
- Rimu.BitStringAddresses.full_mode_maps
- Rimu.BitStringAddresses.hopnextneighbour
- Rimu.BitStringAddresses.near_uniform
- Rimu.BitStringAddresses.num_occupied_modes
- Rimu.BitStringAddresses.occupied_mode_map
- Rimu.BitStringAddresses.occupied_modes
- Rimu.BitStringAddresses.onr
- Rimu.BitStringAddresses.time_reverse
- Rimu.BitStringAddresses.unoccupied_mode_map
- Rimu.BitStringAddresses.unoccupied_modes
- Rimu.DictVectors.append_collections!
- Rimu.DictVectors.append_empty_column!
- Rimu.DictVectors.collect_local!
- Rimu.DictVectors.copy_to_local!
- Rimu.DictVectors.first_column
- Rimu.DictVectors.from_initiator_value
- Rimu.DictVectors.initiator_valtype
- Rimu.DictVectors.is_distributed
- Rimu.DictVectors.local_segments
- Rimu.DictVectors.merge_remote_reductions
- Rimu.DictVectors.move_and_compress!
- Rimu.DictVectors.mpi_exchange_allgather!
- Rimu.DictVectors.mpi_exchange_alltoall!
- Rimu.DictVectors.mpi_recv_any!
- Rimu.DictVectors.mpi_send
- 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.target_segment
- Rimu.DictVectors.to_initiator_value
- Rimu.DictVectors.total_num_segments
- Rimu.DictVectors.walkernumber
- Rimu.DictVectors.walkernumber_and_length
- Rimu.ExactDiagonalization.BasisSetRep
- Rimu.ExactDiagonalization.build_basis
- Rimu.Hamiltonians.HardwallBoundaries
- Rimu.Hamiltonians.LadderBoundaries
- Rimu.Hamiltonians.PeriodicBoundaries
- Rimu.Hamiltonians.check_address_type
- Rimu.Hamiltonians.continuum_dispersion
- Rimu.Hamiltonians.dimension
- Rimu.Hamiltonians.flip_spin_components
- 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.is_invalid_state
- Rimu.Hamiltonians.is_void_state
- Rimu.Hamiltonians.linear_to_state
- Rimu.Hamiltonians.log_abs_oscillator_zero
- Rimu.Hamiltonians.modify_diagonal
- Rimu.Hamiltonians.modify_offdiagonal
- Rimu.Hamiltonians.momentum
- Rimu.Hamiltonians.neighbor_site
- Rimu.Hamiltonians.number_conserving_bose_dimension
- Rimu.Hamiltonians.number_conserving_dimension
- Rimu.Hamiltonians.number_conserving_fermi_dimension
- Rimu.Hamiltonians.one_electron_diagonal
- Rimu.Hamiltonians.rayleigh_quotient
- Rimu.Hamiltonians.shift_lattice
- Rimu.Hamiltonians.shift_lattice_inv
- Rimu.Hamiltonians.two_electron_diagonal
- Rimu.Hamiltonians.unrank_combination
- Rimu.InterfaceTests.test_hamiltonian_interface
- Rimu.InterfaceTests.test_hamiltonian_structure
- Rimu.InterfaceTests.test_observable_interface
- Rimu.InterfaceTests.test_operator_interface
- Rimu.Interfaces.allows_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.dot_from_right
- Rimu.Interfaces.get_offdiagonal
- Rimu.Interfaces.has_adjoint
- Rimu.Interfaces.has_iterable_offdiagonals
- Rimu.Interfaces.has_random_offdiagonal
- Rimu.Interfaces.localpart
- Rimu.Interfaces.num_components
- Rimu.Interfaces.num_modes
- Rimu.Interfaces.num_offdiagonals
- Rimu.Interfaces.num_overlaps
- Rimu.Interfaces.num_particles
- Rimu.Interfaces.num_replicas
- Rimu.Interfaces.num_spectral_states
- Rimu.Interfaces.offdiagonals
- Rimu.Interfaces.operator_column
- Rimu.Interfaces.parent_operator
- Rimu.Interfaces.parent_operator
- Rimu.Interfaces.random_offdiagonal
- Rimu.Interfaces.sort_into_targets!
- Rimu.Interfaces.starting_address
- Rimu.Interfaces.step_stats
- Rimu.Interfaces.storage
- Rimu.Interfaces.sum_mutating!
- Rimu.Interfaces.undo_transform
- Rimu.Interfaces.working_memory
- Rimu.RimuIO.load_df
- Rimu.RimuIO.load_state
- Rimu.RimuIO.save_df
- Rimu.RimuIO.save_state
- Rimu.StatsTools.autocovariance
- Rimu.StatsTools.blocker
- Rimu.StatsTools.blocking_analysis
- Rimu.StatsTools.blocking_analysis_data
- Rimu.StatsTools.blocks_with_m
- Rimu.StatsTools.determine_constant_time_step
- 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.advance!
- Rimu.all_overlaps
- Rimu.clean_and_warn_if_others_present
- Rimu.default_logger
- Rimu.default_starting_vector
- Rimu.delete_and_warn_if_present
- Rimu.finalize_report!
- Rimu.get_metadata
- Rimu.initialise_shift_parameters
- Rimu.is_mpi_root
- Rimu.lomc!
- Rimu.mpi_allprintln
- Rimu.mpi_barrier
- Rimu.mpi_comm
- Rimu.mpi_rank
- Rimu.mpi_rank
- Rimu.mpi_seed!
- Rimu.mpi_size
- Rimu.mpi_size
- Rimu.post_step_action
- Rimu.refine_reporting_strategy
- Rimu.replace_keys
- Rimu.replica_stats
- Rimu.report!
- Rimu.report!
- Rimu.report!
- Rimu.report_after_step!
- Rimu.report_metadata!
- Rimu.reporting_interval
- Rimu.set_up_initial_shift_parameters
- Rimu.single_particle_density
- Rimu.smart_logger
- Rimu.split_keys
- Rimu.state_vectors
- Rimu.update_shift_parameters!
- Rimu.update_time_step
- SparseArrays.sparse
- Statistics.cov
- VectorInterface.scalartype
- Rimu.@mpi_root
- Rimu.BitStringAddresses.@fs_str