API

This page documents the public interface of QILaplace.jl.

MPS Representation

QILaplace.Mps.SignalMPSType
SignalMPS{I<:Index} <: AbstractMPS

A Matrix Product State representing a 1D signal on n qubit sites.

Fields

  • data::Vector{ITensor} — the MPS tensors (one per site).
  • sites::Vector{I} — physical site indices.
  • bonds::Vector{I} — virtual bond indices (length n-1).
  • amplitude::Float64 — the original Euclidean norm of the signal before normalisation.

Use signal_mps to construct a SignalMPS from a dense vector.

source
QILaplace.Mps.ZTMPSType
ZTMPS{I<:Index} <: AbstractMPS

A paired-register MPS for non-unitary transforms (Damping Transform, z-Transform). Each site consists of a PairCore holding a "main" and a "copy" tensor connected by an intra-site bond.

Fields

  • data::Vector{PairCore} — paired tensor cores (one per site).
  • bonds_main::Vector{I} — inter-site bonds (length n-1).
  • bonds_copy::Vector{I} — intra-site bonds (length n).
  • sites_main::Vector{I} — physical site indices on the main register.
  • sites_copy::Vector{I} — physical site indices on the copy register.
  • amplitude::Float64 — the original Euclidean norm of the signal.

Use signal_ztmps to construct a ZTMPS from a dense vector.

source
QILaplace.Mps.PairCoreType
PairCore(Amain::ITensor, Acopy::ITensor, c::Index)

A paired tensor core used inside ZTMPS. Each PairCore holds two ITensors (Amain and Acopy) that share exactly one common index c (the intra-site bond). When a ZTMPS is printed, each site shows its PairCore components.

Fields

  • Amain::ITensor — tensor on the "main" register.
  • Acopy::ITensor — tensor on the "copy" register.
  • c::Index — the shared intra-site bond index.
source

Construction

QILaplace.Signals.generate_signalFunction
generate_signal(n; kind=:sin, dt=nothing, freq=nothing, kwargs...)

Generate a length-2^n real signal of a specified type.

Arguments

  • n::Int: The exponent determining the length of the signal. The signal will have 2^n points.

Keyword Arguments

  • kind::Symbol: The type of signal to generate. Defaults to :sin. Supported kinds include:
    • :sin: A sinusoidal signal.
    • :multi_sin: Multi-frequency sinusoid with deterministic random amplitudes/frequencies.
    • :sin_decay: Exponentially decaying sinusoid (scalar or vector frequency).
    • :multi_sin_exp: Multi-sine with exponential envelopes (deterministic random parameters).
    • :abs_cos_power_p8: |cos(2π*dt*j)|^0.8.
    • :random: A purely random noise signal.
  • dt::Union{Nothing, Float64}: The time step between samples. If nothing, it is automatically computed from freq as dt = 1 / (freq * 2^n). Defaults to nothing.
  • freq::Union{Nothing, Float64, Vector{Float64}}: The frequency or frequencies of the signal. Can be a scalar or a vector of frequencies. Defaults to if not provided.

kwargs (kind-dependent)

  • phase::Float64: The phase offset of the signal in radians. Applicable to :sin and :cos kinds. Defaults to 0.0.
  • decay_rate::Float64: The rate of exponential decay. Only applicable to the :decay kind. A higher value results in faster decay. Defaults to 1.0.
  • noise_level::Float64: The amplitude of random noise added to the signal. Applicable to all kinds. Defaults to 0.0 (no noise).
  • seed::Int: The random seed for reproducibility of noise generation. Applicable when noise_level > 0. Defaults to nothing (random seed).

Returns

  • signal::Vector{Float64}: A real-valued vector of length 2^n representing the generated signal.

Examples

n = 5

# simple sine wave with frequency 1.0
x = generate_signal(n, kind=:sin, freq=1.0) 

# multi-frequency decaying sine wave
x = generate_signal(n, kind=:sin_decay, freq=[1.0, 2.0], decay_rate=[0.1, 0.2])
source
QILaplace.SignalConverters.signal_mpsFunction
signal_mps(x::AbstractVector{<:Number}; method=:svd, kwargs...) -> SignalMPS

Convert a dense signal vector x into a SignalMPS. The vector is zero-padded to the next power of 2 if necessary, normalised, and decomposed into an MPS via SVD (or randomised SVD with method=:rsvd).

The original Euclidean norm is stored in .amplitude.

Keyword Arguments

  • method::Symbol=:svd:svd or :rsvd.
  • cutoff, maxdim — truncation parameters forwarded to the decomposition.

Examples

ψ = signal_mps([1.0, 2.0, 3.0, 4.0])
ψ.amplitude ≈ norm([1,2,3,4])   # true
source
QILaplace.SignalConverters.signal_ztmpsFunction
signal_ztmps(x::AbstractVector{<:Number}; cutoff=1e-10, maxdim=typemax(Int)) -> ZTMPS

Convert a dense signal vector x into a ZTMPS paired-register representation. This duplicates the signal onto a "main" and "copy" register, which is required for non-unitary transforms (Damping Transform, z-Transform).

Keyword Arguments

  • cutoff::Real=1e-10 — SVD truncation cutoff.
  • maxdim::Int — maximum bond dimension.
  • Additional keywords (e.g. method=:rsvd, k, p, q) are forwarded to signal_mps.
source

Extraction & Indexing

QILaplace.Mps.coefficientFunction
coefficient(ψ::SignalMPS, config)
coefficient(ψ::ZTMPS, config)

Return the amplitude ⟨config|ψ⟩ for the given zero-based bit configuration.

config can be:

  • A Vector or Tuple of integers, each in [0, dim(site)-1].
  • A bit string like "1010" or "[1,0,1,0]".
  • A non-negative integer interpreted as an n-bit big-endian pattern.

You may also use direct indexing: ψ[0, 1, 0, 1] is equivalent to coefficient(ψ, (0, 1, 0, 1)).

Examples

ψ = signal_mps([1.0, 0.0, 0.0, 0.0])
coefficient(ψ, [0, 0])   # amplitude of |00⟩
coefficient(ψ, "00")      # same, via bit string
ψ[0, 0]                   # same, via indexing
source
QILaplace.Mps.mps_to_vectorFunction
mps_to_vector(ψ::SignalMPS; reverse::Bool=false)

Extract the dense state vector from a SignalMPS, scaled by the stored amplitude.

Arguments

  • reverse::Bool=false: If false (default), return the vector in MSB-first (normal/natural) bit ordering — matching the original signal ordering from signal_mps. If true, return the vector in the raw column-major (bit-reversed) ordering of the tensor network, which corresponds to the output ordering after a QFT or other transform.

Examples

ψ = signal_mps([1.0, 2.0, 3.0, 4.0])
mps_to_vector(ψ)                  # [1.0, 2.0, 3.0, 4.0] (original signal)
mps_to_vector(ψ; reverse=true)    # bit-reversed ordering
source
mps_to_vector(ψ::ZTMPS; reverse::Bool=false)

Extract the dense state vector from a ZTMPS by first converting to its 2N-site SignalMPS representation, then extracting as a flat vector. Scaled by the stored amplitude.

Arguments

  • reverse::Bool=false: Controls bit ordering (see mps_to_vector(::SignalMPS)).
source
QILaplace.Mps.siteindicesFunction
siteindices(ψ::SignalMPS)
siteindices(ψ::ZTMPS)

Return a named tuple (main=..., copy=...) containing the physical site indices of the MPS. For SignalMPS, copy is an empty vector.

source
QILaplace.Mps.bondindicesFunction
bondindices(ψ::SignalMPS)
bondindices(ψ::ZTMPS)

Return a named tuple (main=..., copy=...) containing the virtual bond indices of the MPS. For SignalMPS, copy is an empty vector.

source

Manipulation

QILaplace.Mps.canonicalize!Function
canonicalize!(ψ, direction; center=nothing, cutoff=1e-12, maxdim=typemax(Int))

Bring the MPS ψ into canonical form (in-place) via QR/LQ sweeps.

Arguments

  • direction::String: "->" (left-canonical) or "<-" (right-canonical).
  • center::Int: orthogonality center site (defaults to N for "->", 1 for "<-").
  • cutoff, maxdim: truncation parameters.

Works for both SignalMPS and ZTMPS.

source
canonicalize!(ψ::ZTMPS, direction::Symbol; center::Union{Nothing,Int}=nothing,
              cutoff::Float64=1e-12, maxdim::Int=typemax(Int))

Bring ZTMPS into canonical form by sweeping through PairCore structures.

Arguments

  • ψ: ZTMPS to canonicalize in-place
  • direction: :right (left-to-right) or :left (right-to-left)
  • center: orthogonality center pair index (defaults to n for :right, 1 for :left)
  • cutoff: truncation threshold for singular values
  • maxdim: maximum bond dimension

Details

For each pair, canonicalizes both Amain and Acopy tensors while maintaining the PairCore structure. Sweeps through pairs in the specified direction.

source
QILaplace.Mps.compress!Function
compress!(ψ; maxdim=typemax(Int), tol=1e-12, sweeps=1)

Truncate bond dimensions of the MPS ψ via alternating SVD sweeps (in-place). After compression, the MPS is re-normalised to unit norm.

Works for both SignalMPS and ZTMPS.

source
QILaplace.Mps.update_site!Function
update_site!(W::SingleSiteMPO, old_site_index::Index, new_site_index::Index)
update_site!(W::PairedSiteMPO, old_site_index::Index, new_site_index::Index)

Replace old_site_index with new_site_index in the MPO W (in-place). Both indices must have the same dimension.

source
update_site!(ψ::SignalMPS, old_site_index::Index, new_site_index::Index)
update_site!(ψ::ZTMPS, old_site_index::Index, new_site_index::Index)

Replace old_site_index with new_site_index in the MPS ψ (in-place). Both indices must have the same dimension.

source
QILaplace.Mps.update_bond!Function
update_bond!(W::SingleSiteMPO, old_bond_index::Index, new_bond_index::Index)
update_bond!(W::PairedSiteMPO, old_bond_index::Index, new_bond_index::Index)

Replace old_bond_index with new_bond_index in the MPO W (in-place). Both indices must have the same dimension.

source
update_bond!(ψ::SignalMPS, old_bond_index::Index, new_bond_index::Index)
update_bond!(ψ::ZTMPS, old_bond_index::Index, new_bond_index::Index)

Replace old_bond_index with new_bond_index in the MPS ψ (in-place). Both indices must have the same dimension.

source
LinearAlgebra.normFunction
norm(ψ::SignalMPS)
norm(ψ::ZTMPS)

Compute the Euclidean norm √⟨ψ|ψ⟩ of the MPS state by full tensor contraction. This extends LinearAlgebra.norm.

source

MPO Representation

QILaplace.Mpo.SingleSiteMPOType
SingleSiteMPO{I<:Index} <: AbstractMPO

A Matrix Product Operator acting on single-register MPS (e.g. QFT). Each tensor has physical indices s and s' plus virtual bond indices.

Fields

  • data::Vector{ITensor} — the MPO tensors.
  • sites::Vector{I} — physical site indices.
  • bonds::Vector{I} — virtual bond indices (length n-1).
source
QILaplace.Mpo.PairedSiteMPOType
PairedSiteMPO{I<:Index} <: AbstractMPO

A Matrix Product Operator acting on paired-register MPS (ZTMPS). Used for the Damping Transform and z-Transform. The operator's 2n tensors alternate between main and copy sites.

Fields

  • data::Vector{ITensor} — the MPO tensors (length 2n).
  • sites_main, sites_copy — physical site indices on each register.
  • bonds_main, bonds_copy — inter-site and intra-site bond indices.
source

Transformers

QILaplace.QFTTransform.build_qft_mpoFunction
build_qft_mpo(n, sites; cutoff=1e-14, maxdim=1000) -> SingleSiteMPO
build_qft_mpo(ψ::SignalMPS; cutoff=1e-14, maxdim=1000)

Build the Quantum Fourier Transform MPO for n qubits. The MPO is composed from controlled-Hadamard-phase gates assembled via the zip-up / zip-down compression algorithm. The cutoff tells the error threshold of the compressed MPO.

When called with a SignalMPS, the site indices are taken from the MPS itself.

Examples

ψ = signal_mps(rand(16))
W_qft = build_qft_mpo(ψ)
ψ_freq = apply(W_qft, ψ)   # frequency-domain MPS
source
QILaplace.DTTransform.build_dt_mpoFunction
build_dt_mpo(n, ωr, sites_main, sites_copy; cutoff=1e-14, maxdim=1000) -> PairedSiteMPO
build_dt_mpo(ψ::ZTMPS, ωr; cutoff=1e-14, maxdim=1000)

Build the Damping Transform MPO for n qubits with damping parameter ωr. This applies an exponential envelope exp(-ωr * k) in the MPS representation, corresponding to the damping part of the discrete Laplace transform. The cutoff tells the error threshold of the compressed MPO.

The MPO is composed of two parts built via zip-to-combine / zip-to-compress:

  1. Control-damping blocks on main sites (k = 1 → n).
  2. Control-damping-copy blocks on copy sites (k = n → 2).

When called with a ZTMPS, site indices are taken from the MPS itself.

Examples

ψ_z = signal_ztmps(x)
W_dt = build_dt_mpo(ψ_z, ωr)
ψ_out = apply(W_dt, ψ_z)
source
QILaplace.ZTTransformer.build_zt_mpoFunction
build_zt_mpo(n, ωr, sites_main, sites_copy; cutoff=1e-14, maxdim=1000) -> PairedSiteMPO
build_zt_mpo(ψ::ZTMPS, ωr; cutoff=1e-14, maxdim=1000)

Build the discrete Laplace transform (z-Transform) MPO for n qubits at damping parameter ωr. The cutoff tells the error threshold of the compressed MPO.

Operationally, this builds the DT (build_dt_mpo) and the paired ZTMPS QFT (control_Hphase_ztmps_mpo on all 2n sites), then fuses them with apply(W_dt, W_qft) and compresses — same end operator as the legacy interleaved construction, without per-stage DT×QFT growth.

When called with a ZTMPS, site indices are taken from the MPS itself.

Examples

ψ_z = signal_ztmps(x)
W_zt = build_zt_mpo(ψ_z, ωr)
ψ_out = apply(W_zt, ψ_z)
source

Linear Algebra Utilities

ITensors.applyFunction
apply(W, ψ; kwargs...) -> ψ_out
W * ψ

Contract an MPO W with an MPS ψ, returning a new MPS.

Supported dispatch combinations:

  • SingleSiteMPO × SignalMPS — standard MPO-MPS contraction.
  • PairedSiteMPO × ZTMPS — paired-register contraction (Damping / z-Transform).
  • SingleSiteMPO × SingleSiteMPO — MPO composition.
  • PairedSiteMPO × PairedSiteMPO — paired MPO composition.

The * operator is an alias for apply.

source
QILaplace.RSVD.rsvdFunction
rsvd(A::ITensor, Linds...; k=20, p=10, q=0, random_seed=1234, bondtag="Link,rsvd",
     verbose=false, cutoff=1e-15, maxdim=k, mindim=1) -> (U, S, V)

Randomised SVD for an ITensor, returning factors A ≈ U * S * V where S is diagonal. Uses native ITensor qr/svd internally to avoid dense matrix conversions.

Arguments

  • A — the ITensor to decompose.
  • Linds — indices that should end up on the left factor U.

Keyword Arguments

  • k::Int=20 — target rank (number of kept singular values).
  • p::Int=10 — oversampling; the projection uses k + p random vectors.
  • q::Int=0 — number of power iterations (improves accuracy on slowly decaying spectra).
  • random_seed::Int — seed for the random projection.
  • cutoff::Float64=1e-15 — singular value truncation threshold.
  • maxdim::Int=k — hard cap on bond dimension.
  • mindim::Int=1 — minimum bond dimension to keep.
  • bondtag — tag string applied to the SVD bond index.
  • verbose::Bool — print power-iteration progress.

When to use

Use :rsvd (via signal_mps(x; method=:rsvd, k=...)) when the signal is very large and you want a fast, low-rank approximation instead of a full SVD.

source