API

Rimu

Rimu.RimuModule
Rimu

Random integrators for many-body quantum systems

Welcome to Rimu version 0.14.2. Read the documentation online.

source
Rimu.AllOverlapsType
AllOverlaps(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}.

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

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

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

Note that adjusting the keyword maxlength in ProjectorMonteCarloProblem is advised as the default may not be appropriate.

See ShiftStrategy, ProjectorMonteCarloProblem.

source
Rimu.DoubleLogSumUpdateType
DoubleLogSumUpdate(; target_walkers = 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, ProjectorMonteCarloProblem.

source
Rimu.DoubleLogUpdateType
DoubleLogUpdate(; target_walkers = 1_000, ζ = 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, ProjectorMonteCarloProblem.

source
Rimu.FCIQMCType
FCIQMC(; kwargs...) <: PMCAlgorithm

Algorithm 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.

source
Rimu.FirstOrderTransitionOperatorType
FirstOrderTransitionOperator(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.

source
Rimu.GramSchmidtType
GramSchmidt(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.

source
Rimu.MultiScalarType
MultiScalar

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

Example

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

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

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

julia> s, p
(0, 2)

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

source
Rimu.PMCSimulationType
PMCSimulation

Holds the state and the results of a projector quantum Monte Carlo (PMC) simulation. Is returned by init(::ProjectorMonteCarloProblem) and solved with solve!(::PMCSimulation).

Obtain the results of a simulation sm as a DataFrame with DataFrame(sm).

Fields

  • problem::ProjectorMonteCarloProblem: The problem that was 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!.

source
Rimu.PostStepStrategyType
PostStepStrategy

Subtypes 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).

source
Rimu.ProjectedEnergyType
ProjectedEnergy(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.

source
Rimu.ProjectorMonteCarloProblemType
ProjectorMonteCarloProblem(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 StochasticStyle of 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 true to force PDVec for the starting vector, false for 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 > 1 to find excited states.
  • spectral_strategy = GramSchmidt(n_spectral): The SpectralStrategy used for orthogonalizing spectral states.

Example

julia> hamiltonian = HubbardReal1D(BoseFS(1,2,3));

julia> problem = ProjectorMonteCarloProblem(hamiltonian; target_walkers = 500, last_step = 100);

julia> simulation = solve(problem);

julia> simulation.success[]
true

julia> size(DataFrame(simulation))
(100, 9)

Further keyword arguments:

  • starting_step = 1: Starting step of the simulation.
  • wall_time = Inf: Maximum time allowed for the simulation.
  • simulation_plan = SimulationPlan(; starting_step, last_step, wall_time): Defines the duration of the simulation. Takes precedence over last_step and 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 FCIQMC is implemented.
  • shift: Initial shift value or collection of shift values. Determined by default from the Hamiltonian and the starting vectors.
  • initial_shift_parameters: Initial shift parameters or collection of initial shift parameters. Overrides shift if 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 ProjectorMonteCarloProblem twice 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_at is not provided.

See also init, solve.

source
Rimu.ReplicaStateType
ReplicaState <: 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.

source
Rimu.ReplicaStrategyType
ReplicaStrategy{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:

source
Rimu.ReportDFAndInfoType
ReportDFAndInfo(; reporting_interval=1, info_interval=100, io=stdout, writeinfo=false) <: ReportingStrategy

The 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.

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

ReportingStrategy 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 DataFrame of 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 DataFrame is returned.
  • io = stdout: The IO to print messages to. Set to devnull if you don't want to see messages printed out.
  • compress = :zstd: compression algorithm to use. Can be :zstd, :lz4 or nothing.

See also load_df, save_df, ReportDFAndInfo, and ProjectorMonteCarloProblem.

source
Rimu.RunTillLastStepType
RunTillLastStep(step::Int = 0 # number of current/starting timestep
             laststep::Int = 100 # number of final timestep
             shiftMode::Bool = false # whether to adjust shift
             shift = 0.0 # starting/current value of shift
             dτ::Float64 = 0.01 # current value of time step
) <: FciqmcRunStrategy

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

Warning

The use of this strategy is deprecated. Pass the relevant arguments directly to ProjectorMonteCarloProblem or to lomc!() instead.

source
Rimu.ShiftStrategyType
ShiftStrategy

Abstract 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:

Extended help

Internally To implement a custom strategy, define a new subtype of Rimu.ShiftStrategy and implement methods for:

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

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

source
Rimu.SingleParticleDensityType
SingleParticleDensity(; 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

source
Rimu.SingleStateType
SingleState(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.

source
Rimu.SpectralStrategyType
SpectralStrategy{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.
source
Rimu.TimerType
Timer <: PostStepStrategy

Record current time after every step. See Base.Libc.time for information on what time is recorded.

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

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

source
CommonSolve.solve!Method
solve!(sm::PMCSimulation; kwargs...)::PMCSimulation

Solve a Rimu.PMCSimulation until the last step is completed or the wall time limit is reached.

To continue a previously completed simulation, set a new last_step or wall_time using the keyword arguments. Optionally, changes can be made to the replica_strategy, the post_step_strategy, or the reporting_strategy.

Optional keyword arguments:

  • last_step = nothing: Set the last step to a new value and continue the simulation.
  • wall_time = nothing: Set the allowed wall time to a new value and continue the simulation.
  • reset_time = false: Reset the elapsed_time counter and continue the simulation.
  • empty_report = false: Empty the report before continuing the simulation.
  • replica_strategy = nothing: Change the replica strategy. Requires the number of replicas to match the number of replicas in the 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.

source
Rimu.advance!Method
advance!(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.

See also solve!, step!.

source
Rimu.all_overlapsMethod
all_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}.

source
Rimu.clean_and_warn_if_others_presentMethod
clean_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.

source
Rimu.default_loggerMethod
default_logger(args...)

Reset the global_logger to Logging.ConsoleLogger. Undoes the effect of smart_logger. Arguments are passed on to Logging.ConsoleLogger.

source
Rimu.default_starting_vectorMethod
default_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.

source
Rimu.delete_and_warn_if_presentMethod
delete_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.

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

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

source
Rimu.lomc!Method
lomc!(ham::AbstractHamiltonian, [v]; kwargs...) -> df, state
lomc!(state::ReplicaState, [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 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 v and shift.
  • style = IsStochasticInteger() - set StochasticStyle for default v; unused if v is specified.
  • initiator = NonInitiator() - set InitiatorRule for default v; unused if v is specified.
  • threading - default is to use multithreading and MPI if multiple threads are available. Set to true to force PDVec for the starting vector, false for serial computation; unused if v is 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_step is accepted as an alias for post_step_strategy.)
  • replica_strategy::ReplicaStrategy = NoStats(1) - run several synchronised simulations, see ReplicaStrategy. (Deprecated: replica is accepted as an alias for replica_strategy.)
  • reporting_strategy::ReportingStrategy = ReportDFAndInfo() - how and when to report results, see ReportingStrategy. (Deprecated: r_strat is 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 DataFrame with all statistics being reported.
  • state: a ReplicaState that 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 AbstractHamiltonian argument, a DataFrame can 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.

Warning

The use of this lomc! is deprecated. Use ProjectorMonteCarloProblem and solve instead.

source
Rimu.mpi_allprintlnMethod
mpi_allprintln(args...)

Print a message to stdout from each rank separately, in order. MPI synchronizing.

source
Rimu.mpi_barrierFunction
mpi_barrier(comm = mpi_comm())

The MPI barrier with optional argument. MPI syncronizing.

source
Rimu.mpi_seed!Function
mpi_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.

source
Rimu.refine_reporting_strategyMethod
refine_reporting_strategy(reporting_strategy::ReportingStrategy) -> reporting_strategy

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

source
Rimu.replace_keysMethod
replace_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.

source
Rimu.replica_statsFunction
replica_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.

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

Report keys and values to report, which will be converted to a DataFrame before 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...).

source
Rimu.report!Method
report!(report::Report, df::DataFrame)

Convert the DataFrame df to a Report. This function does not copy the data.

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

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

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

source
Rimu.report_metadata!Method
report_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.

source
Rimu.reporting_intervalMethod
reporting_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.

source
Rimu.set_up_initial_shift_parametersMethod
set_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.

source
Rimu.single_particle_densityMethod
single_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

source
Rimu.smart_loggerMethod
smart_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().

source
Rimu.split_keysMethod
split_keys(nt::NamedTuple, selected_keys...) -> extracted, rest

Find 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))
source
Rimu.state_vectorsMethod
state_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.

source
Rimu.update_time_stepMethod
update_time_step(s<:TimeStepStrategy, time_step, tnorm) -> new_time_step

Update the time step according to the strategy s.

source
Rimu.@mpi_rootMacro
@mpi_root expr

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

Example:

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

Reexported Submodules

ExactDiagonalization

See Exact Diagonalization

Interfaces

See Module Interfaces

StochasticStyles

See Module StochasticStyles

Hamiltonians

See Module Hamiltonians

BitStringAddresses

See Module BitStringAddresses

DictVectors

See Module DictVectors

StatsTools

See Module StatsTools

Index