API reference

ReducedBasis.AffineDecompositionType

Represents an affine decomposition

\[O(\bm{\mu}) = \sum_{r=1}^R \alpha_r(\bm{\mu})\, O_r\]

where terms<:AbstractArray carries the $O_r$ matrices (or more generally linear maps) and the coefficients implements the $\alpha_r(\bm{\mu})$ coefficient functions.

In the case where the $\alpha_r$ are constant coefficients, i.e. that don't depend on the parameter vector, coefficients can also just be an AbstractArray of the same size as terms.

Note that $r = (r_1, \dots, r_d)$ generally is a multi-index and as such terms can be a $d$-dimensional array. Correspondingly, coefficients maps parameter points $\bm{\mu}$ to an size(terms) array.

source
ReducedBasis.ApproxMPOType

Carries an ITensors.MPO matrix-product operator and possible truncation keyword arguments.

This enables a simple mpo * mps syntax while allowing for proper truncation throughout the basis assembly. Includes the exact operator sum in opsum to be able to produce efficient sums of MPOs when constructing AffineDecomposition sums explicitly.

source
ReducedBasis.ApproxMPOMethod
ApproxMPO(mpo::MPO, opsum; <keyword arguments>)

Construct an ApproxMPO with default truncation settings.

Arguments

  • mpo::MPO
  • opsum::Sum{Scaled{ComplexF64,Prod{Op}}}
  • cutoff::Float64=1e-9: relative cutoff for singular values.
  • maxdim::Int=1000: maximal bond dimension.
  • mindim::Int=1: minimal bond dimension.
source
ReducedBasis.DMRGType

Solver type for the density matrix renormalization group (DMRG) as implemented in ITensors.dmrg.

Fields

  • n_states::Int=1: see FullDiagonalization.
  • tol_degeneracy::Float64=0.0: see FullDiagonalization.
  • sweeps::Sweeps=default_sweeps(; cutoff_max=1e-9, bonddim_max=1000): set DMRG sweep settings via ITensors.Sweeps.
  • observer::Function=() -> DMRGObserver(; energy_tol=1e-9): set DMRG exit conditions. At each solve a new ITensors.AbstractObserver object is created.
  • verbose::Bool=false: if true, prints info about DMRG solve.
source
ReducedBasis.EigenDecompositionType

Extension type for orthogonalization and compression using eigenvalue decomposition of the basis overlap matrix. See also extend.

Fields

  • cutoff::Float64=1e-6: cutoff for minimal eigenvalue accuracy.
source
ReducedBasis.FullDiagonalizationType

Solver type for full diagonalization using LinearAlgebra.eigen.

Fields

  • n_target::Int=1: the number of the targeted eigenvalue. If n_target=1, the degeneracy of the ground state is automatically determined up to tol_degeneracy. If tol_degeneracy=0, it determines the number of returned vectors per solve (also includes excited states).
  • tol_degeneracy::Float64=0.0: tolerance for distinguishing two eigenvalues. If abs(λ₁ - λ₂) < tol_degeneracy, the eigenvalues are added to the same degenerate subspace.
source
ReducedBasis.GreedyType

Greedy reduced basis assembling strategy.

Fields

  • estimator::ErrorEstimate: error estimate used in greedy condition. See also estimate_error
  • tol::Float64=1e-3: tolerance for error estimate, below which the assembly is terminated.
  • n_truth_max::Int=64: maximal number of truth solves to be taken up in the basis.
  • Ψ_init::Function=rb_guess: returns initial guess for truth solver. See also interpolate.
  • verbose::Bool=true: print information during assembly if true.
  • exit_checks::Bool=true: if false, no exit checks will be performed and the assembly will run until tol or n_truth_max are reached.
source
ReducedBasis.HamiltonianCacheType

Convenience type storing a Hamiltonian, its applications to vectors and its compressions.

Fields

  • H::AffineDecomposition
  • HΨ::Vector{V}
  • ΨHΨ::Vector{Matrix{T}}
  • ΨHHΨ::Matrix{Matrix{T}}
  • h::AffineDecomposition
  • h²::AffineDecomposition
source
ReducedBasis.InfoCollectorMethod
(collector::InfoCollector)(info)

Push iteration information into InfoCollector and return info object, containing the InfoCollector itself.

source
ReducedBasis.InfoCollectorMethod
InfoCollector(fields::Symbol...)

Construct InfoCollector from fields that are contained in the info iteration state object. Possible fields to select from are:

  • iteration: number of iteration at which the information was obtained.
  • err_grid: error estimate on all parameter points of the training grid.
  • λ_grid: RB energies on all training grid points.
  • err_max: maximal error estimate on the grid.
  • μ: parameter point at which truth solve has been performed.
  • solver_info: output of the solving method, which excludes eigenvalues and vectors.
  • basis: RBasis at the current iteration.
  • extend_info: info that is specific to the chosen extension procedure, not including

the extended RBasis.

  • condnum: condition number of the $B^\dagger B$.
  • h: current reduced Hamiltonian.
  • h_cache: HamiltonianCache at the current iteration.
source
ReducedBasis.LOBPCGType

Solver type for the Locally Optimal Block Preconditioned Conjugate Gradient (LOBPCG). Currently uses the DFTK lobpcg_hyper implementation.

Fields

  • n_target::Int=1: see FullDiagonalization.
  • tol_degeneracy::Float64=0.0: see FullDiagonalization.
  • tol::Float64=1e-9: tolerance for residual norms.
  • maxiter::Int=300: maximal number of LOBPCG iterations.
  • n_ep_extra::Int=4: number of extra eigenpairs that are kept to improve convergence.
  • shift::Float64=-100: eigenvalue shift.
  • verbose::Bool=false: if true, print convergence messages.
  • dense_fallback::Bool=true: if false, also non-converged states will be accepted. Otherwise LinearAlgebra.eigen is used for non-converged iterations.
  • maxdiagonal::Int=400
source
ReducedBasis.PODType

Proper orthogonal decomposition assembly strategy.

Fields

  • n_vectors::Int: number of retained singular vectors of the snapshot matrix in the returned basis.
  • verbose::Bool=true: shows the truth solve progress if true.
source
ReducedBasis.QRCompressType

Extension type for QR orthonormalization and compression. See extend for details.

Fields

  • tol::Float64: tolerance for compressing insignificant basis snapshots.
source
ReducedBasis.RBasisType

Central type containing snapshots and associated objects that make a reduced basis.

The snapshot vectors are contained in snapshots::AbstractVector{V} where the snapshots are of type V. Here a "snapshot" refers to all vectors that were obtained from a solve at a specific parameter point. For greedy-generated bases, each element in snapshots has a parameter point in parameter::Vector{P}, i.e. for snapshots $\bm{\Psi}(\bm{\mu_i}) = (\Psi_1(\bm{\mu_i}),\dots,\Psi_m(\bm{\mu_i}))$ of multiplicity $m$ the parameter point $\bm{\mu}$ is contained $m$ times. However, for general basis-assembly strategies this one-to-one correspondence might not be true, e.g. for POD where snapshots contains singular vectors based on all truth solves.

Treated as a matrix, the reduced basis $B = \Upsilon V$ is made up of the snapshot vectors as column vectors in $\Upsilon$ and vector coefficients $V$. The latter are stored in vectors. In the simple case of $B = \Upsilon$, one sets vectors=I.

Since the matrix $B^\dagger B = V^\dagger \Upsilon^\dagger \Upsilon V$ is frequently needed, both the snapshot_overlaps $\Upsilon^\dagger \Upsilon$ and the metric $B^\dagger B$ are stored with generic floating-point type T.

source
ReducedBasis.ResidualType

Estimator type for the residual $\mathrm{Res}(\bm{\mu}) = \lVert H(\bm{\mu}) B \varphi(\bm{\mu}) - \lambda B \varphi(\bm{\mu}) \rVert$.

source
Base.:*Method
*(o::ApproxMPO, mps::MPS)

Apply o.mpo to mps using the truncation arguments contained in o.

source
Base.truncateMethod
Base.truncate(ad_raw::AffineDecomposition, basis_trunc::RBasis)
Base.truncate(ad::AffineDecomposition,
              basis_trunc::RBasis{V,T,P,<:UniformScaling}) where {V,T<:Number,P}

Truncate an AffineDecomposition.

For general RBasis, the "raw" AffineDecomposition has to be provided, where ad_raw.terms are the matrix elements $\Psi_i^\dagger O_r \Psi_j$. The returned decomposition then contains the fully compressed transformed terms, i.e. $B^\dagger O_r B$ using the truncated basis.

In case a basis is provided that has trivial vectors=I, the truncation is performed naively on the fully compressed terms.

source
Base.truncateMethod
Base.truncate(hc::HamiltonianCache, basis_trunc::RBasis)

Truncate a [HamiltonianCache] according to an already truncated basis.

source
Base.truncateMethod
Base.truncate(basis::RBasis, n_truth::Int)

Truncate the RBasis to n_truth snapshots.

Note that n_truth does not amount to the dimension of the truncated basis, but the number of truth solves included in the basis, which can feature degenerate snapshots.

source
ReducedBasis.assembleMethod
assemble(H::AffineDecomposition, grid, pod::POD, solver_truth)

Assemble basis using POD.

Only ED solvers such as FullDiagonalization and LOBPCG are supported. The generated RBasis will contain pod.n_vectors singular vectors in snapshots and all grid points in parameters. This means that parameters and snapshots generally have different lengths.

source
ReducedBasis.assembleMethod
assemble(H, grid, greedy, solver_truth, compressalg; <keyword arguments>)
assemble(info, H, grid, greedy, solver_truth, compressalg; <keyword arguments>)

Assemble an RBasis using the greedy strategy and any truth solving method.

If info is not provided, start a greedy assembly from scratch.

Arguments

  • info::NamedTuple: iteration state from a previous simulation that will be resumed.
  • H::AffineDecomposition: Hamiltonian for which a reduced basis is assembled.
  • grid: parameter grid that defines the parameter domain.
  • greedy::Greedy: greedy strategy containing assembly parameters. See also Greedy.
  • solver_truth: solving method for obtaining ground state snapshots.
  • compressalg: compression method for orthogonalization, etc. See also extend.
  • solver_online=FullDiagonalization(solver_truth): solving method that is used for the RB generalized eigenvalue problem.
  • callback=print_callback: callback function that operates on the iteration state during assembly. It is possible to chain multiple callback functions using .
  • μ_start=grid[1]: parameter point of the starting iteration, if no info is provided.
source
ReducedBasis.compressMethod
compress(M::AbstractMatrix, snapshots::AbstractVector{<:AbstractVector})

Compress term using matrix multiplications.

source
ReducedBasis.compressMethod
compress(ad::AffineDecomposition{<:Matrix,<:Function},
         basis::RBasis; symmetric_terms=false)

Perform compression for an AffineDecomposition with terms with two indices (double-sum observables), including an option to exploit the possible symmetry of terms $O_{r,r'} = O_{r',r}$, such that only the necessary compressions are computed.

source
ReducedBasis.default_sweepsMethod
default_sweeps(; cutoff_max=1e-9, bonddim_max=1000, iter_max=100)

Return default ITensors.Sweeps object for DMRG solves, containing decreasing noise and increasing maximal bond dimension ramps.

source
ReducedBasis.extendMethod
extend(basis::RBasis, new_snapshot::AbstractVector, μ, ed::EigenDecomposition)

Extend the reduced basis by orthonormalizing and compressing via eigenvalue decomposition.

The overlap matrix $S$ in basis.snapshot_overlaps is eigenvalue decomposed $S = U^\dagger \Lambda U$ and orthonormalized by computing the vector coefficients $V = U \Lambda^{-1/2}$. Modes with an relative squared eigenvalue error smaller than ed.cutoff are dropped.

source
ReducedBasis.extendMethod
extend(basis::RBasis, new_snapshot, μ, ::NoCompress)

Extend the reduced basis by one snapshot without any orthogonalization or compression procedure.

source
ReducedBasis.extendMethod
extend(basis::RBasis, new_snapshot::AbstractVector, μ, qrcomp::QRCompress)

Extend using QR orthonormalization and compression.

The orthonormalization is performed by QR decomposing the orthogonal projection $\bm{\Psi}(\bm{\mu}_{n+1}) - B_n^\dagger [B_n^\dagger B_n]^{-1} B_n \bm{\Psi}(\bm{\mu}_n)$ and appending $Q$ to snapshots. Modes that have an $R$ column maximum falling below the qrcomp.tol tolerance are dropped.

source
ReducedBasis.extend_overlapsMethod
extend_overlaps(old_overlaps::Matrix, old_snapshots::Vector, new_snapshot::Vector)

Extend an overlap matrix by one snapshot, where only the necessary dot products are computed.

source
ReducedBasis.in_boundsMethod
in_bounds(μ, grid::RegularGrid)

Check whether a given parameter point μ is in the convex hull of the grid.

source
ReducedBasis.interpolateMethod
interpolate(basis::RBasis, h::AffineDecomposition, μ, solver_online)

Compute Hilbert-space-dimensional ground state vector at parameter point μ from the reduced basis using $\bm{\Phi}(\bm{\mu}) = B \bm{\varphi}(\bm{\mu})$.

source
ReducedBasis.interpolateMethod
interpolate(basis::RBasis{MPS}, h::AffineDecomposition, μ, dm::DMRG, solver_online)

Compute ground state MPS at μ from the reduced basis by MPS addition in $| \Phi(\bm{\mu}) \rangle = \sum_{k=1}^{\dim B} [V \varphi(\bm{\mu_k})]_k\, | \Psi(\bm{\mu_k}) \rangle$.

source
ReducedBasis.overlap_matrixMethod
overlap_matrix(v1::Vector, v2::Vector)
overlap_matrix(f, v1::Vector, v2::Vector)

Compute the overlap matrix of two sets of vector-like objects v1 and v2.

The computed matrix elements are the dot products dot(v1[i], v2[j]). Correspondingly, the elements of v1 and v2 must support a LinearAlgebra.dot method. In the case where v1 = v2, the Gram matrix is computed.

Optionally, a function f can be provided which is applied to each element of v2.

source
ReducedBasis.rb_guessMethod
rb_guess(info, offline_args; kwargs...)

Provide the reduced basis interpolated ground state as an initial guess for the truth solver.

source
ReducedBasis.reconstructMethod
reconstruct(mps::MPS)

Explicitly compute Hilbert-space-dimensional vector by reconstructing all MPS coefficients.

Memory restrictions

The number of MPS coefficients grows exponentially with system size, such that the explicit reconstruction is only possible for small systems.

source
ReducedBasis.shiftMethod
shift(grid::RegularGrid{D,N}, μ_shift; stay_in_bounds=false) where {D,N}

Shift a regular grid by a shift vector μ_shift. If stay_in_bounds=true, the shifted grid will stay in the convex hull of the unshifted grid. Note that the shift vector elements cannot be larger than the grid.ranges steps.

source
ReducedBasis.solveMethod
solve(H::AffineDecomposition, μ, Ψ₀::Matrix, lobpcg::LOBPCG)
solve(H::AffineDecomposition, μ, ::Nothing, lobpcg::LOBPCG)

Solve using LOBPCG. If nothing is provided as an initial guess, an orthogonal random matrix will be used with lobpcg.n_target + lobpcg.n_ep_extra column vectors.

source
ReducedBasis.solveMethod
solve(H::AffineDecomposition, μ, Ψ₀::Vector{MPS}, dm::DMRG)
solve(H::AffineDecomposition, μ, Ψ₀::MPS, dm::DMRG)
solve(H::AffineDecomposition, μ, ::Nothing, dm::DMRG)

Solve using DMRG.

The length of the Ψ₀ vector determines the number of targeted states, given that dm.n_states > 1 and dm.tol_degeneracy > 0. When nothing is provided as an initial guess, dm.n_states random MPS are used.

source
ReducedBasis.solveMethod
solve(h::AffineDecomposition, b::Matrix, μ, fd::FullDiagonalization)

Solve the generalized eigenvalue problem $h(\bm{\mu})\, \varphi(\bm{\mu}) = \lambda(\bm{\mu})\, b\, \varphi(\bm{\mu})$ at parameter point μ using FullDiagonalization.

source