Skip to main content

dmrg

Overloads

NameDescription
dmrg(const ComplexMPS initial_state, const ComplexMPO hamiltonian, const Map options) -> ListCompute the ground state of a quantum many-body system using the Density Matrix Renormalization Group (DMRG) algorithm.
dmrg(const RealMPS initial_state, const RealMPO hamiltonian, const Map options) -> ListCompute the ground state of a quantum many-body system using the Density Matrix Renormalization Group (DMRG) algorithm.

dmrg(const ComplexMPS initial_state, const ComplexMPO hamiltonian, const Map options) -> List

Compute the ground state of a quantum many-body system using the Density Matrix Renormalization Group (DMRG) algorithm.

DMRG is a variational method that optimizes a Matrix Product State (MPS) representation to minimize the energy expectation value with respect to a given Hamiltonian (MPO). This implementation uses a two-site update algorithm with adaptive bond dimension truncation.

The algorithm supports multiple optimization schedules and configurable logging for monitoring convergence.

Parameters

  • initial_state: Initial MPS guess for the ground state. Must have the same number of sites as the Hamiltonian. The algorithm will optimize this state to minimize energy.
  • hamiltonian: Hamiltonian operator in MPO form representing the quantum system.
  • options: Map of algorithm parameters.
    • schedule ("static", "logarithmic", "polynomial"): Optimization schedule type controlling bond dimension growth. default: "logarithmic". Formulas (χ\chi = current max bond dimension): • "static": χ\chi remains constant (fixed max bond dimension). • "logarithmic": χnew=χ+log2(χ+1)×slope\chi_{\text{new}} = \chi + \lfloor \log_2(\chi + 1) \times \text{slope} \rfloor (slope 1.0\geq 1.0) • "polynomial": χnew=χ+slope×χ(degree1)/degree\chi_{\text{new}} = \chi + \lfloor \text{slope} \times \chi^{(\text{degree}-1)/\text{degree}} \rfloor (slope 0.0\geq 0.0)
    • max_iterations (positive integer): Maximum number of DMRG sweeps. default: 100.
    • energy_accuracy (positive real number): Energy convergence tolerance. Stops when |E_new - E_old| < tolerance. default: 1e-04.
    • truncation_accuracy (positive real number): SVD truncation tolerance for bond dimension reduction. default: 1e-10.
    • max_bond_dimension (positive integer): Maximum allowed bond dimension. default: 18446744073709551615.
    • degree (positive integer): Polynomial degree (only for polynomial schedule). default: 1.
    • slope (positive real number): Schedule slope parameter affecting convergence rate. default: 2.
    • verbose ("none", "sweep", "update", "all"): Logging level for monitoring convergence. default: "none".

Returns

List containing two elements:

  • [0]: Ground state energy (real number)
  • [1]: Optimized ground state MPS

Example

// Compute ground state of 1D transverse-field Ising model
var L = 4
var H = Ising_1d(L, complex(0.5)) // transverse field strength = 0.5, using complex numbers
var psi0 = make_random_mps(L, 2, 2, as_complex)

// Configure DMRG with static schedule and minimal logging
var options = [
"schedule": "static",
"max_iterations": 10,
"max_bond_dimension": 4,
"energy_accuracy": 1.0e-8,
"truncation_accuracy": 1.0e-8,
"verbose": "sweep"
]

var result = dmrg(psi0, H, options)
var energy = result[0]
var ground_state = result[1]
print("Ground state energy: " + string(energy))

dmrg(const RealMPS initial_state, const RealMPO hamiltonian, const Map options) -> List

Compute the ground state of a quantum many-body system using the Density Matrix Renormalization Group (DMRG) algorithm.

DMRG is a variational method that optimizes a Matrix Product State (MPS) representation to minimize the energy expectation value with respect to a given Hamiltonian (MPO). This implementation uses a two-site update algorithm with adaptive bond dimension truncation.

The algorithm supports multiple optimization schedules and configurable logging for monitoring convergence.

Parameters

  • initial_state: Initial MPS guess for the ground state. Must have the same number of sites as the Hamiltonian. The algorithm will optimize this state to minimize energy.
  • hamiltonian: Hamiltonian operator in MPO form representing the quantum system.
  • options: Map of algorithm parameters.
    • schedule ("static", "logarithmic", "polynomial"): Optimization schedule type controlling bond dimension growth. default: "logarithmic". Formulas (χ\chi = current max bond dimension): • "static": χ\chi remains constant (fixed max bond dimension). • "logarithmic": χnew=χ+log2(χ+1)×slope\chi_{\text{new}} = \chi + \lfloor \log_2(\chi + 1) \times \text{slope} \rfloor (slope 1.0\geq 1.0) • "polynomial": χnew=χ+slope×χ(degree1)/degree\chi_{\text{new}} = \chi + \lfloor \text{slope} \times \chi^{(\text{degree}-1)/\text{degree}} \rfloor (slope 0.0\geq 0.0)
    • max_iterations (positive integer): Maximum number of DMRG sweeps. default: 100.
    • energy_accuracy (positive real number): Energy convergence tolerance. Stops when |E_new - E_old| < tolerance. default: 1e-04.
    • truncation_accuracy (positive real number): SVD truncation tolerance for bond dimension reduction. default: 1e-10.
    • max_bond_dimension (positive integer): Maximum allowed bond dimension. default: 18446744073709551615.
    • degree (positive integer): Polynomial degree (only for polynomial schedule). default: 1.
    • slope (positive real number): Schedule slope parameter affecting convergence rate. default: 2.
    • verbose ("none", "sweep", "update", "all"): Logging level for monitoring convergence. default: "none".

Returns

List containing two elements:

  • [0]: Ground state energy (real number)
  • [1]: Optimized ground state MPS

Example

// Compute ground state of 1D transverse-field Ising model
var L = 4
var H = Ising_1d(L, real(0.5)) // transverse field strength = 0.5, using real numbers
var psi0 = make_random_mps(L, 2, 2, as_real)

// Configure DMRG with static schedule and minimal logging
var options = [
"schedule": "static",
"max_iterations": 10,
"max_bond_dimension": 4,
"energy_accuracy": 1.0e-8,
"truncation_accuracy": 1.0e-8,
"verbose": "sweep"
]

var result = dmrg(psi0, H, options)
var energy = result[0]
var ground_state = result[1]
print("Ground state energy: " + string(energy))