Module DictVectors

Rimu.Interfaces.AbstractDVecType
AbstractDVec{K,V}

Abstract data type for vector-like data structures with sparse storage. While conceptually AbstractDVecs represent elements of a vector space over a scalar type V, they are indexed by an arbitrary type K (could be non-integers) similar to dictionaries. They support the interface from VectorInterface.jl and are designed to work well for quantum Monte Carlo with ProjectorMonteCarloProblem and for matrix-free linear algebra with KrylovKit.

Concrete implementations are available as PDVec, DVec, and InitiatorDVec.

AbstractDVecs have a StochasticStyle which selects the spawning algorithm in FCIQMC. Looking up an element that is not stored in the AbstractDVec should return a zero, and setting a value to zero should remove it from the vector. To iterate over an AbstractDVec, use keys, pairs, or values. When possible, use reduction functions such as sum or mapreduce.

Interface

The interface is similar to the AbstractDict interface, but with the changed behaviour as noted above. Implement what would be needed for the AbstractDict interface (pairs, keys, values, setindex!, getindex, delete!, length, empty, empty!) and, in addition:

A default implementation for the VectorInterface.jl interface is provided through the above functions.

See also DictVectors, Interfaces.

source

Concrete implementations

Rimu.DictVectors.DVecType
DVec(args...; kwargs...) <: AbstractDVec{K,V}

Dictionary-based vector-like data structure for storing key <: K and value <: V pairs for use with ProjectorMonteCarloProblem and ExactDiagonalizationProblem. DVec supports the interface from VectorInterface.jl but stores only non-zero values.

DVec is fast but does not support Initiators (see InitiatorDVec) or parallel and distributed vector operations (see PDVec).

See also: AbstractDVec, InitiatorDVec, PDVec.

Constructors

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

  • DVec(args...[; style, capacity]): args... can be key => value Pairs or an iterator over Pairs or Tuples.

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

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

Keyword arguments

  • style::StochasticStyle determines the mode of stochastic operations. See StochasticStyle. The default style is selected based on the DVec's valtype (see default_style). If a style is given and the valtype does not match the style's eltype, the values are converted to an appropriate type.
  • capacity: The capacity argument is optional and sets the initial size of the DVec via Base.sizehint!.

Examples

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

julia> dv = DVec(:a => 2, :b => 3; style=IsDeterministic())
DVec{Symbol,Float64} with 2 entries, style = IsDeterministic{Float64}()
  :a => 2.0
  :b => 3.0
source
Rimu.DictVectors.InitiatorDVecType
InitiatorDVec(args...; kwargs...) <: AbstractDVec{K,V}

Dictionary-based vector-like data structure for storing key <: K and value <: V pairs for use with ProjectorMonteCarloProblem and KrylovKit.jl. See AbstractDVec. Functionally identical to DVec, but contains InitiatorValues internally in order to facilitate initiator methods. Initiator methods for controlling the Monte Carlo sign problem were first introduced in J. Chem. Phys. 132, 041103 (2010). How the initiators are handled is controlled by specifying an InitiatorRule with the initiator keyword argument (see below). Note that values need to be comparable to real numbers with isless.

See also: AbstractDVec, DVec, PDVec.

Constructors

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

  • InitiatorDVec(args...[; style, initiator, capacity]): args... can be key => value Pairs or an iterator over Pairs or Tuples.

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

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

Keyword arguments

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

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

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

source
Rimu.DictVectors.PDVecType
PDVec{K,V}(; kwargs...) <: AbstractDVec{K,V}
PDVec(iter; kwargs...)
PDVec(pairs...; kwargs...)

Dictionary-based vector-like data structure for storing key <: K and value <: V pairs for use with ProjectorMonteCarloProblem and ExactDiagonalizationProblem. Only pairs with non-zero values are stored. PDVec supports various linear algebra operations such as norm and dot, as well as the full interface defined in VectorInterface.

The P in PDVec stands for parallel. PDVecs perform reductions such as mapreduce, foreach, and various linear algebra operations in a threaded manner on the available Julia threads. If MPI is available, the storage is distributed over the MPI ranks and reductions operate in parallel on the full distributed data structure. As such it is not recommended to iterate over pairs, keys, or values directly unless explicitly performing them on the localpart of the vector.

See also: AbstractDVec, DVec, InitiatorDVec.

Keyword arguments

Example

julia> address = FermiFS2C((1,1,0,0), (0,0,1,1));

julia> op = HubbardMom1D(address; t=4/π^2, u=4);

julia> pv = PDVec(address => 1.0)
1-element PDVec: style = IsDeterministic{Float64}()
  fs"|↑↑↓↓⟩" => 1.0

julia> pv = op * pv
7-element PDVec: style = IsDeterministic{Float64}()
  fs"|↑↑↓↓⟩" => 4.0
  fs"|↑↓↑↓⟩" => 1.0
  fs"|⋅⇅⇅⋅⟩" => 1.0
  fs"|↓↑↑↓⟩" => -1.0
  fs"|↑↓↓↑⟩" => -1.0
  fs"|⇅⋅⋅⇅⟩" => 1.0
  fs"|↓↑↓↑⟩" => 1.0

julia> scale!(pv, -1); pv
7-element PDVec: style = IsDeterministic{Float64}()
  fs"|↑↑↓↓⟩" => -4.0
  fs"|↑↓↑↓⟩" => -1.0
  fs"|⋅⇅⇅⋅⟩" => -1.0
  fs"|↓↑↑↓⟩" => 1.0
  fs"|↑↓↓↑⟩" => 1.0
  fs"|⇅⋅⋅⇅⟩" => -1.0
  fs"|↓↑↓↑⟩" => -1.0

julia> dest = similar(pv)
0-element PDVec: style = IsDeterministic{Float64}()

julia> map!(x -> x + 2, dest, values(pv))
7-element PDVec: style = IsDeterministic{Float64}()
  fs"|↑↑↓↓⟩" => -2.0
  fs"|↑↓↑↓⟩" => 1.0
  fs"|⋅⇅⇅⋅⟩" => 1.0
  fs"|↓↑↑↓⟩" => 3.0
  fs"|↑↓↓↑⟩" => 3.0
  fs"|⇅⋅⋅⇅⟩" => 1.0
  fs"|↓↑↓↑⟩" => 1.0

julia> sum(values(pv))
-6.0

julia> dot(dest, pv)
10.0

julia> dot(dest, op, pv)
44.0

Extended Help

Segmentation

The vector is split into Threads.nthreads() subdictionaries called segments. Which dictionary a key-value pair is mapped to is determined by the hash of the key. The purpose of this segmentation is to allow parallel processing - functions such as mapreduce, add! or dot (full list below) process each subdictionary on a separate thread.

See also PDWorkingMemory.

MPI

When MPI is active, all parallel reductions are automatically reduced across MPI ranks with a call to MPI.Allreduce!.

In a distributed setting, PDVec does not support iteration without first making it explicit the iteration is only to be performed on the local segments of the vector. This is done with localpart. In general, even when not using MPI, it is best practice to use localpart when explicit iteration is required.

Use with KrylovKit

PDVec is compatible with eigsolve from KrylovKit.jl. When used, the diagonalisation is performed in a threaded and distributed manner. Using multiple MPI ranks with this method does not distribute the memory load effectively, but does result in significant speedups.

Example

julia> using KrylovKit

julia> address = BoseFS((0,0,5,0,0));

julia> op = HubbardMom1D(address; u=6.0);

julia> pv = PDVec(address => 1.0);

julia> results = eigsolve(op, pv, 4, :SR);

julia> results[1][1:4]
4-element Vector{Float64}:
 -3.431115689232223
  1.182174860261238
  3.737775375308286
  6.996390417443125

Parallel functionality

The following functions are threaded and MPI-compatible:

  • From Base: mapreduce and derivatives (sum, prod, reduce...), all, any,map! (on values only), +, -, *
  • From LinearAlgebra: rmul!, lmul!, mul!, axpy!, axpby!, dot, norm, normalize, normalize!
  • The full interface defined in VectorInterface.jl
source

Interface functions

Rimu.Interfaces.deposit!Function
deposit!(w::InitiatorDVec, add, val, p_add=>p_val)

Add val into w at address add as an AbstractInitiatorValue.

source
deposit!(w::AbstractDVec, add, val, parent::Pair)

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

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

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

source
OrderedCollections.freezeFunction
freeze(dv)

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

source
Rimu.Interfaces.apply_operator!Function
apply_operator!(working_memory, target, source, operator, boost=1, compress=Val(true)) ->
    stat_names, stats, working_memory, target

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

\[v^{(n + 1)} = \hat{T} v^{(n)} ,\]

where $\hat{T}$ is the operator, $v^{(n+1)}$ is the target and $v^{(n)}$ is the source. The working_memory can be used as temporary storage.

The boost argument is passed to apply_column! and increases the number of spawns performed. For the operator to be applied without compressing the vector after, set compress to Val(false).

Whether the operation is performed in a stochastic, semistochastic, or determistic way is controlled by the trait StochasticStyle(target). See StochasticStyle.

Returns the step stats generated by the StochasticStyle, the working memory and the target vector. target and working_memory may be mutated and/or swapped.

source
Rimu.Interfaces.sort_into_targets!Function
sort_into_targets!(target, source, stats) -> target, source, agg_stats

Aggregate coefficients from source to target and from stats to agg_stats according to thread- or MPI-level parallelism.

Returns the new target and source, as the sorting process may involve swapping them.

source
Base.mapreduceFunction
mapreduce(f, op, keys(::PDVec); [init])
mapreduce(f, op, values(::PDVec); [init])
mapreduce(f, op, pairs(::PDVec); [init])

Perform a parallel reduction operation on PDVecs. MPI-compatible. Is used in the definition of various functions from Base such as reduce, sum, prod, etc.

init, if provided, must be a neutral element for op.

source
Rimu.Interfaces.sum_mutating!Function
sum_mutating!(accumulator, [f! = add!], keys(::PDVec); [init])
sum_mutating!(accumulator, [f! = add!], values(::PDVec); [init])
sum_mutating!(accumulator, [f! = add!], pairs(::PDVec); [init])

Perform a parallel sum on PDVecs for vector-valued results while minimizing allocations. The result of the sum will be added to accumulator and stored in accumulator. MPI-compatible. If f! is provided, it must accept two arguments, the first being the accumulator and the second the element of the iterator. Otherwise,add! is used.

If provided, init must be a neutral element for + and of the same type as accumulator.

See also mapreduce.

source
sum_mutating!(accumulator, [f! = add!], iterator)

Add the sum of elements in iterator to accumulator, storing the result in accumulator. If f! is provided, it must accept two arguments, the first being the accumulator and the second the element of the iterator. Otherwise, add! is used.

See also mapreduce.

source

Supported operations

AbstractDVecs generally support most operations that are defined on Vectors and Dicts. This includes the interface from VectorInterface.jl, and many functions from the LinearAlgebra standard library.

A significant difference between AbstractDVecs, Vectors, and Dicts, is that iteration on them is disabled by default. Iteration must be explicitly performed on keys, values, or pairs, however, it is highly recommended you use mapreduce, reduce, or similar functions when performing reductions, as that will make the operations compatible with MPI.

In addition, Rimu defines the following function.

Rimu.DictVectors.walkernumberFunction
walkernumber(v)

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

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

source

Projectors

Initiator rules

Rimu.DictVectors.InitiatorType
Initiator(threshold = 1.0) <: InitiatorRule

Initiator rule to be passed to PDVec or InitiatorDVec. An initiator is a configuration add with a coefficient with magnitude abs(v[add]) > threshold. The threshold can be passed as a keyword argument. Rules:

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

See InitiatorRule.

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

Initiator rule to be passed to PDVec or InitiatorDVec. An initiator is a configuration add with a coefficient with magnitude abs(v[add]) > threshold. The threshold can be passed as a keyword argument. Rules:

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

See InitiatorRule.

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

Initiator rule to be passed to PDVec or InitiatorDVec. An initiator is a configuration add with a coefficient with magnitude abs(v[add]) > threshold. The threshold can be passed as a keyword argument. Rules:

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

See InitiatorRule.

source

PDVec internals

Working memory

Rimu.DictVectors.PDWorkingMemoryType
PDWorkingMemory(t::PDVec)

The working memory that handles threading and MPI distribution for operations that involve operators, such as FCIQMC propagation, operator-vector multiplication and three-way dot products with PDVecs.

The working memory is structured as a two-dimensional array of segments, which themselves are Dicts (see PDVec). The number of rows in this array is equal to the number of segments across all MPI ranks (covering the entire address space), while the number of columns corresponds to the number of segments in the current MPI rank (i.e. column corresponds to the part of the address space that is local to the current rank).

The purpose of this organisation is to allow spawning in parallel without using locks or atomic operations. The spawning is performed by applying the following sequence of operations:

  • perform_spawns!: each segment in the PDVec is multiplied by the operator independently, with the results being stored in a column of the working memory.
  • collect_local!: the rows of the working memory are summed to the first column.
  • synchronize_remote!: the segments corresponding to other MPI ranks are distributed and transferred to the first column.
  • move_and_compress!: the results are stochastically compressed and moved to the result PDVec

When used with three-argument dot products, a full copy of the left-hand side vector is materialized in the first column of the working memory on all ranks.

source

Communicators

Rimu.DictVectors.AllToAllType
AllToAll{K,V}(; mpi_comm, n_segments, report) <: Communicator

Communicator that uses collective communication using MPI.Alltoall[v]!.

Keyword arguments

  • mpi_comm=MPI.COMM_WORLD: the MPI communicator to use.
  • n_segments=Threads.nthreads(): the number of segments per rank to use. Should match the PDVec the communicator is used with.
  • report=false: if set to true, report MPI communication times during a projector Monte Carlo run.

See also: Communicator.

source
Rimu.DictVectors.CommunicatorType
abstract type Communicator

Communicators are used to handle MPI communication when using PDVecs. Currently, three implementations are provided, NotDistributed, AllToAll and PointToPoint. The communicator is picked automatically according to the number of MPI ranks available.

When implementing a communicator, use local_segments and remote_segments.

Interface

Optional interface

See also: PDVec, PDWorkingMemory.

source
Rimu.DictVectors.NestedSegmentedBufferType
NestedSegmentedBuffer{T}(nrows) <: AbstractMatrix{AbstractVector{T}}

Matrix of vectors stored in a single buffer with collective MPI communication support. The number of rows in the matrix is fixed to nrows.

Used in the AllToAll communication strategy, where each column corresponds to an MPI rank and each row corresponds to a segment in the PDVec.

Supported operations

See also: SegmentedBuffer.

source
Rimu.DictVectors.PointToPointType
PointToPoint{K,V}(; mpi_comm, report) <: Communicator

MPI Communicator that uses circular communication using MPI.Isend and MPI.Recv!.

Keyword arguments

  • mpi_comm=MPI.COMM_WORLD: the MPI communicator to use.
  • report=false: if set to true, report MPI communication times during a projector Monte Carlo run.
source
Rimu.DictVectors.SegmentedBufferType
SegmentedBuffer{T}() <: AbstractVector{AbstractVector{T}}

Behaves like a vector of vectors, but is stored in a single buffer. It can be sent/received over MPI keeping its structure intact. Used in the PointToPoint communication strategy.

Supported operations

See also: NestedSegmentedBuffer.

source
Rimu.DictVectors.mpi_recv_any!Method
mpi_recv_any!(buf::SegmentedBuffer, comm::MPI_Comm) -> Int

Find a source that is ready to send a buffer and receive from it. Return the rank ID of the sender.

source
Rimu.DictVectors.replace_collections!Method
replace_collections!(buf::SegmentedBuffer, iters)

Insert collections in iters into a SegmentedBuffer.

julia> using Rimu.DictVectors: SegmentedBuffer

julia> buf = SegmentedBuffer{Int}()
0-element SegmentedBuffer{Int64}

julia> Rimu.DictVectors.replace_collections!(buf, [[1,2,3], [4,5]])
2-element SegmentedBuffer{Int64}:
 [1, 2, 3]
 [4, 5]

julia> Rimu.DictVectors.replace_collections!(buf, [[1], [2,3], [4]])
3-element SegmentedBuffer{Int64}:
 [1]
 [2, 3]
 [4]
source
Rimu.DictVectors.target_segmentMethod
target_segment(c::Communicator, k, num_segments) -> target, is_local

This function is used to determine where in the PDVec a key should be stored. If the key is local (stored on the same MPI rank), return its segment index and true. If the key is non-local, return any value and false.

See also: PDVec, Communicator.

source

Index