ScarFinder
ScarFinder in MPSToolkit.jl is a deliberately explicit workflow:
- evolve the state
- project or truncate it
- optionally match a target energy
- optionally refine along a short projected trajectory
That split is important. ScarFinder is not a hidden backend with its own internal state machine. It is a small, composable loop built from public pieces like evolve!, project!, match_energy!, and selector scoring.
When To Use It
ScarFinder is useful when you want to search for low-entanglement or approximately recurrent trajectories without baking that search logic directly into the evolution engine itself. In practice, that often means:
- evolve under a physical or effective Hamiltonian for a short time
- project back into a chosen variational manifold
- repeat until the trajectory stabilizes
- rank or refine the resulting orbit with a simple selector
Recommended Step-Count Rule
ScarFinder now treats an effective evolution step count of 1 as a bad main-loop setting.
- If
scarfinder_step!,scarfinder!, ortrajectory_refine!receive:LocalGateEvolutionwithnstep == 1DMTGateEvolutionwithnstep == 1TDVPEvolutionwith effective steps1vianstepsornsweeps
- then ScarFinder emits a warning and internally uses
10for that ScarFinder call.
This rule is local to ScarFinder only. The global TEBD, DMT, and TDVP constructors still keep their own defaults. If you already know you are building an evolution object for ScarFinder, it is best to set nstep=10 or nsteps=10 explicitly and avoid the warning.
Internal correction moves inside match_energy! still use one-step updates on purpose. Those are narrow post-step corrections, not the main ScarFinder trajectory evolution.
Minimal Workflow
using MPSToolkit
truncation = BondDimTruncation(16; cutoff=1e-10)
target = EnergyTarget(0.0; operator=hamiltonian_mpo, tol=1e-8, alpha=0.05, maxstep=8)
selector = EntropySelector(; bond=3)
scarfinder!(
psi,
evolution,
truncation,
10;
target_energy=target,
refine=true,
selector=selector,
refine_steps=6,
)The important point is that projection is explicit. scarfinder_step! does not hide a bespoke evolution backend; it composes the same public building blocks that you can also call yourself.
PXP Background
The PXP model is the standard kinetically constrained spin-chain model for Rydberg-blockaded atom arrays, and it is one of the best-known settings for quantum many-body scar dynamics. In the ScarFinder paper, this projection loop is used on the PXP model to uncover a trajectory with nearly perfect revivals in the thermodynamic limit without assuming prior knowledge of the scar tower.
For an open chain, the Hamiltonian used in the example is
\[H_{\mathrm{PXP}} = \sum_{j=2}^{L-1} 2 P^{\downarrow}_{j-1} X_j P^{\downarrow}_{j+1},\]
where P^\downarrow = |\downarrow\rangle\langle\downarrow| enforces the blockade constraint on neighboring sites and X_j flips the center spin.
Conceptually, the PXP use case looks like:
- choose a low-entanglement initial product-state ansatz
- evolve it for a short time under the constrained Hamiltonian
- explicitly project back to the chosen variational manifold
- iterate until the trajectory converges toward a stable revival orbit
That makes PXP a good mental model for what ScarFinder is for: not just finding low-entanglement states, but isolating coherent trajectories embedded inside otherwise thermalizing many-body dynamics.
Example
- examples/scarfinder/pxp_scarfinder.ipynb A
L = 32open-chain PXP notebook using 3-site TEBD gates, explicit projection, energy targeting, and selector-based diagnostics.
Core PXP Setup
The notebook uses the public API directly. A minimal ScarFinder-oriented setup is:
using MPSToolkit
using ITensors
using ITensorMPS
using LinearAlgebra
function pxp_local_hamiltonian()
projector_dn = ComplexF64[0 0; 0 1]
pauli_x = ComplexF64[0 1; 1 0]
return kron(projector_dn, kron(pauli_x, projector_dn))
end
function pxp_schedule(nsites::Int)
starts = Int[]
for offset in 1:3
append!(starts, offset:3:(nsites - 2))
end
return starts
end
function pxp_mpo(sites)
opsum = OpSum()
for j in 2:(length(sites) - 1)
opsum += 2.0, "ProjDn", j - 1, "Sx", j, "ProjDn", j + 1
end
return MPO(opsum, sites)
end
nsites = 32
sites = siteinds("S=1/2", nsites)
schedule = pxp_schedule(nsites)
local_hamiltonian = pxp_local_hamiltonian()
hamiltonian_mpo = pxp_mpo(sites)
scar_evolution = tebd_evolution_from_hamiltonians(
fill(local_hamiltonian, length(schedule)),
0.1;
schedule=schedule,
nstep=10,
maxdim=64,
cutoff=1e-10,
)
truncation = BondDimTruncation(1; cutoff=1e-10)
energy_target = EnergyTarget(0.0; operator=hamiltonian_mpo, tol=1e-8, alpha=0.1, maxstep=4)
z2 = MPS(sites, n -> isodd(n) ? "Up" : "Dn")
psi = deepcopy(z2)
for _ in 1:200
scarfinder_step!(psi, scar_evolution, truncation; target_energy=energy_target)
endIf you are using TDVP instead of dense-gate TEBD, the same recommendation applies: prefer an explicit nsteps=10 configuration when the evolution object is intended for ScarFinder.
Selectors And Targeting
BondDimTruncationcontrols the explicit projection budget.EnergyTargetadds a post-step correction loop toward a desired expectation value.EntropySelectorandFidelitySelectorprovide small scoring rules for trajectory refinement.SelectionContextcarries optional external information, such as the reference state needed byFidelitySelector.
This makes the search logic easy to inspect and easy to modify. If you want to replace the projector, selector, or energy-correction step, you can do that without rewriting the entire workflow.
API
MPSToolkit.BondDimTruncation — Type
BondDimTruncationTruncation settings used by ScarFinder projection routines.
Fields
maxdim: Maximum allowed bond dimension after projection.cutoff: Singular-value cutoff used during truncation.
MPSToolkit.EnergyTarget — Type
EnergyTargetParameters controlling post-evolution energy correction toward a target expectation value.
Fields
target: Desired energy density or expectation value.operator: Dense operator orMPOused to measure and correct the energy.tol: Early-stop tolerance on the residual energy mismatch.alpha: Proportional step-size parameter for the correction update.maxstep: Maximum number of correction iterations per ScarFinder step.
MPSToolkit.SelectionContext — Type
SelectionContextAdditional reference data passed to selector scoring during ScarFinder refinement.
Fields
reference_state: Optional state used by selectors such asFidelitySelector.
MPSToolkit.EntropySelector — Type
EntropySelectorSelector configuration that scores states using bond entanglement entropy.
Fields
bond: Bond whose entropy should be minimized.nothingmeans "use the backend default".
MPSToolkit.FidelitySelector — Type
FidelitySelectorSelector configuration that scores states using fidelity against a reference state.
Notes
- The required reference state is supplied separately through
SelectionContext.
MPSToolkit.project! — Function
project!(psi, truncation)Backend-dispatched in-place projection or truncation entry point.
Arguments
psi: Mutable state to project back into the chosen variational manifold.truncation: Backend-specific truncation or projection settings.
Returns
- The same
psiobject after in-place mutation.
Notes
- ScarFinder uses this after each evolution step, which keeps the projection logic explicit instead of hiding it inside the evolution backend.
project!(psi::MPS, trunc)Truncate a finite MPS in place using standard ITensor/ITensorMPS compression.
Arguments
psi: FiniteMPSto truncate.trunc:BondDimTruncationobject providingmaxdimandcutoff.
Returns
- The mutated
psi.
Notes
- This is the default ScarFinder projection rule for finite-state workflows.
- Compression is delegated to
truncate!, so orthogonality handling follows the underlying ITensor implementation.
MPSToolkit.match_energy! — Function
match_energy!(psi, evolution, truncation, target)Apply a post-evolution correction loop that nudges the state toward a target energy.
Arguments
psi: State to mutate in place.evolution: Main ScarFinder evolution object.truncation: Projection settings reused after each correction substep.target:EnergyTargetconfiguration.
Returns
- The mutated
psi.
Notes
- Dense operators use the local-gate correction path;
MPOs use a one-step TDVP correction path.
MPSToolkit.trajectory_refine! — Function
trajectory_refine!(psi, evolution, truncation, selector; kwargs...)Public ScarFinder refinement entry point.
Notes
- This wrapper applies ScarFinder step-count normalization before dispatching to the internal refinement routine.
MPSToolkit.scarfinder_step! — Function
scarfinder_step!(psi, evolution, truncation; target_energy=nothing, selector=nothing, kwargs...)Run one public ScarFinder step.
Arguments
psi: State to mutate in place.evolution: Evolution object for the main ScarFinder update.truncation: Projection settings applied after evolution.
Keyword Arguments
target_energy: OptionalEnergyTargetused for post-step correction.selector: Accepted for interface symmetry withscarfinder!. It is not used by a single-step update.kwargs...: Reserved for forward compatibility with higher-level workflow wrappers.
Returns
- The mutated
psi.
Notes
- If the effective step count of
evolutionis1, ScarFinder emits a warning and uses10internally for this call.
MPSToolkit.scarfinder! — Function
scarfinder!(psi, evolution, truncation, niter; refine=false, selector=nothing, kwargs...)Run a multi-step ScarFinder loop.
Arguments
psi: State to mutate in place.evolution: Evolution object for the main ScarFinder updates.truncation: Projection settings applied after each step.niter: Number of projected evolution steps to perform.
Keyword Arguments
refine: Iftrue, run selector-based trajectory refinement after the main loop.selector: Selector configuration used during refinement.target_energy: OptionalEnergyTargetreused on every step.kwargs...: Additional refinement keywords forwarded totrajectory_refine!, such asrefine_stepsorselector_context.
Returns
- The mutated
psi.
Notes
- The ScarFinder single-step normalization rule is applied once up front, so repeated loop iterations do not emit repeated warnings.
References
- Jie Ren, Andrew Hallam, Lei Ying, and Zlatko Papic, ScarFinder: A detector of optimal scar trajectories in quantum many-body dynamics