Time Evolution using Renormalizer#
Overview#
In this notebook we will simulate the charge transfer between two molecules using the Marcus model
with transfer integral \(V=-0.1\), dimensionless coupling constant \(g=1\), vibration frequency \(\omega=0.5\).
We will first show how to use Renormalizer to simulation the time evolution at \(\Delta G = 0\). Then, we will show that, by decreasing the reaction Gibbs free energy change \(\Delta G\), the reaction rate will first increase and then decrease, as predicted by the Marcus theory.
Preparation: Setting Up Logger#
[1]:
from renormalizer.utils.log import package_logger as logger
2025-08-18 06:51:00,084[INFO] Use NumPy as backend
2025-08-18 06:51:00,085[INFO] numpy random seed is 9012
2025-08-18 06:51:00,085[INFO] random seed is 1092
2025-08-18 06:51:00,095[INFO] Git Commit Hash: c080e702a6ee0a055fc180ecc6ad885b0368b1c1
2025-08-18 06:51:00,096[INFO] use 64 bits
[2]:
logger.debug("logger output")
2025-08-18 06:51:00,134[DEBUG] logger output
[3]:
from renormalizer.utils.log import set_stream_level, INFO
[4]:
# filter logger output
set_stream_level(INFO)
[5]:
logger.debug("This message will not be shown")
[6]:
logger.info("Logger output")
2025-08-18 06:51:00,155[INFO] Logger output
Define the Model and Initial State#
[7]:
from renormalizer import Op, BasisMultiElectron, BasisSHO, Model
import numpy as np
[8]:
v = -0.1
g = 1
omega = 0.5
nbas = 16
[9]:
def get_model(delta_g):
ham_terms = v * Op(r"a^\dagger a", ["e0", "e1"]) + v * Op(r"a^\dagger a", ["e1", "e0"])
ham_terms += delta_g * Op(r"a^\dagger a", "e1")
for i in range(2):
ham_terms += omega * Op(r"b^\dagger b", f"v{i}")
ham_terms += g * omega * Op(r"a^\dagger a", f"e{i}") * Op(r"b^\dagger+b", f"v{i}")
basis = [BasisMultiElectron(["e0", "e1"], [0, 0]), BasisSHO("v0", omega, nbas), BasisSHO("v1", omega, nbas)]
return Model(basis, ham_terms)
[10]:
# using a relaxed initial state
def get_init_condition():
basis = BasisSHO(0, omega, nbas)
state = np.linalg.eigh(basis.op_mat(r"b^\dagger b") + g * basis.op_mat(r"b^\dagger+b"))[1][:, 0]
return {"v0": state}
init_condition = get_init_condition()
Time Evolution with the Default Configuration#
Next, we run the simulation using the evolve
method in the Mps
class. At this phase, we keep \(\Delta G = 0\).
[11]:
delta_g = 0
model = get_model(delta_g)
[12]:
from renormalizer import Mps, Mpo
[13]:
# Hamiltonian MPO
mpo = Mpo(model)
# The occupation of e0
n_op = Mpo(model, Op(r"a^\dagger a", "e0"))
[14]:
# initialize the MPS
mps = Mps.hartree_product_state(model, condition=init_condition)
# time evolution step
dt = 0.2
# record the electronic occupation
n_list = []
for i_step in range(50):
n = mps.expectation(n_op)
logger.info(f"Step {i_step}. Time {i_step * dt:.2f}. $n$ {n}")
# perform time evolution. Note that the evolution is not in-place.
mps = mps.evolve(mpo, dt)
n_list.append(n)
2025-08-18 06:51:00,205[INFO] Step 0. Time 0.00. $n$ 0.9999999999999998
2025-08-18 06:51:00,219[INFO] Step 1. Time 0.20. $n$ 0.9996030537865594
2025-08-18 06:51:00,231[INFO] Step 2. Time 0.40. $n$ 0.998431396656162
2025-08-18 06:51:00,245[INFO] Step 3. Time 0.60. $n$ 0.9965532876450579
2025-08-18 06:51:00,260[INFO] Step 4. Time 0.80. $n$ 0.9940731362639054
2025-08-18 06:51:00,276[INFO] Step 5. Time 1.00. $n$ 0.9911161868842491
2025-08-18 06:51:00,292[INFO] Step 6. Time 1.20. $n$ 0.9878132450238745
2025-08-18 06:51:00,308[INFO] Step 7. Time 1.40. $n$ 0.9842873078294337
2025-08-18 06:51:00,324[INFO] Step 8. Time 1.60. $n$ 0.9806441927926955
2025-08-18 06:51:00,340[INFO] Step 9. Time 1.80. $n$ 0.976967746624379
2025-08-18 06:51:00,356[INFO] Step 10. Time 2.00. $n$ 0.9733182226560808
2025-08-18 06:51:00,372[INFO] Step 11. Time 2.20. $n$ 0.9697373081106376
2025-08-18 06:51:00,388[INFO] Step 12. Time 2.40. $n$ 0.9662497016753406
2025-08-18 06:51:00,404[INFO] Step 13. Time 2.60. $n$ 0.962867160945706
2025-08-18 06:51:00,420[INFO] Step 14. Time 2.80. $n$ 0.9595919026133962
2025-08-18 06:51:00,436[INFO] Step 15. Time 3.00. $n$ 0.9564197521679629
2025-08-18 06:51:00,452[INFO] Step 16. Time 3.20. $n$ 0.9533423734970323
2025-08-18 06:51:00,469[INFO] Step 17. Time 3.40. $n$ 0.9503492089387037
2025-08-18 06:51:00,485[INFO] Step 18. Time 3.60. $n$ 0.9474287391824951
2025-08-18 06:51:00,501[INFO] Step 19. Time 3.80. $n$ 0.9445692798814913
2025-08-18 06:51:00,518[INFO] Step 20. Time 4.00. $n$ 0.9417594196434944
2025-08-18 06:51:00,534[INFO] Step 21. Time 4.20. $n$ 0.9389883140168839
2025-08-18 06:51:00,551[INFO] Step 22. Time 4.40. $n$ 0.9362456942950406
2025-08-18 06:51:00,567[INFO] Step 23. Time 4.60. $n$ 0.9335221984937244
2025-08-18 06:51:00,584[INFO] Step 24. Time 4.80. $n$ 0.9308095732941701
2025-08-18 06:51:00,601[INFO] Step 25. Time 5.00. $n$ 0.9281008849250423
2025-08-18 06:51:00,617[INFO] Step 26. Time 5.20. $n$ 0.9253906107085343
2025-08-18 06:51:00,634[INFO] Step 27. Time 5.40. $n$ 0.9226738471654728
2025-08-18 06:51:00,650[INFO] Step 28. Time 5.60. $n$ 0.9199478805518423
2025-08-18 06:51:00,667[INFO] Step 29. Time 5.80. $n$ 0.9172094618331406
2025-08-18 06:51:00,684[INFO] Step 30. Time 6.00. $n$ 0.9144551950424458
2025-08-18 06:51:00,700[INFO] Step 31. Time 6.20. $n$ 0.9116816732036944
2025-08-18 06:51:00,717[INFO] Step 32. Time 6.40. $n$ 0.9088859026726781
2025-08-18 06:51:00,733[INFO] Step 33. Time 6.60. $n$ 0.9060659656334908
2025-08-18 06:51:00,750[INFO] Step 34. Time 6.80. $n$ 0.9032210121837198
2025-08-18 06:51:00,767[INFO] Step 35. Time 7.00. $n$ 0.9003515623906606
2025-08-18 06:51:00,784[INFO] Step 36. Time 7.20. $n$ 0.8974590790192766
2025-08-18 06:51:00,801[INFO] Step 37. Time 7.40. $n$ 0.8945453849965326
2025-08-18 06:51:00,818[INFO] Step 38. Time 7.60. $n$ 0.8916120979060871
2025-08-18 06:51:00,835[INFO] Step 39. Time 7.80. $n$ 0.8886602437062094
2025-08-18 06:51:00,852[INFO] Step 40. Time 8.00. $n$ 0.8856902107426247
2025-08-18 06:51:00,868[INFO] Step 41. Time 8.20. $n$ 0.8827020259679245
2025-08-18 06:51:00,885[INFO] Step 42. Time 8.40. $n$ 0.8796958667524066
2025-08-18 06:51:00,901[INFO] Step 43. Time 8.60. $n$ 0.8766726737398559
2025-08-18 06:51:00,917[INFO] Step 44. Time 8.80. $n$ 0.8736347299848458
2025-08-18 06:51:00,935[INFO] Step 45. Time 9.00. $n$ 0.8705861034366976
2025-08-18 06:51:00,951[INFO] Step 46. Time 9.20. $n$ 0.8675328737556515
2025-08-18 06:51:00,968[INFO] Step 47. Time 9.40. $n$ 0.8644830817098144
2025-08-18 06:51:00,984[INFO] Step 48. Time 9.60. $n$ 0.8614463370235249
2025-08-18 06:51:01,001[INFO] Step 49. Time 9.80. $n$ 0.8584330043902285
[15]:
# plotting the occupation
from matplotlib import pyplot as plt
plt.style.use("mm.mplstyle")
plt.plot(np.arange(len(n_list)) * dt, n_list, label="Default")
plt.xlabel("$t$")
plt.ylabel(r"Occupation $\langle a^\dagger_0 a_0 \rangle$")
plt.legend()
plt.show()

In this example, the time evolution is performed with default time evolution configuration and MPS compression configuration.
Renormalizer by default uses the RK4 “propagation and compression” method for time evolution. The advantage of the method is that it is easy to understand and setup.
Please refer to our recent review for a discussion of the different time evolution schemes.
[16]:
mps.evolve_config.method
[16]:
<EvolveMethod.prop_and_compress: 'P&C'>
Regarding MPS compression/truncation configuration, Renormalizer by default employs a truncation scheme based on the the singular value threshold.
[17]:
mps.compress_config.criteria, mps.compress_config.threshold
[17]:
(<CompressCriteria.threshold: 'threshold'>, 0.001)
Inspection of the dimension of the final mps after time evolution shows that the compression is efficient.
[18]:
mps.bond_dims
[18]:
[1, 2, 4, 1]
Configuring Time Evolution#
In order to configure the time evolution, you should update the evolve_config
and compress_config
attributes of an MPS instance. These attributes are instances of the EvolveConfig
class and the CompressConfig
class. Please see the API referece for full options of the configuration classes.
[19]:
from renormalizer.utils.configs import CompressConfig, CompressCriteria
[20]:
# update the compresssion configuration. Adopt a fixed bond dimension of 8
mps.compress_config = CompressConfig(CompressCriteria.fixed, max_bonddim=8)
Next, we perform one more step of the time evolution, and the bond dimension in the middle of the MPS is increased from 4 to 8, according to the updated compression configuration.
[21]:
new_mps = mps.evolve(mpo, dt)
new_mps.bond_dims
[21]:
[1, 2, 8, 1]
Recent studies have shown that methods based on time dependent variantional principle (TDVP) show higher accuracy and efficiency. In our production runs, we usually employ one-site TDVP with projector splitting (TDVP-PS) for time evolution. TDVP-PS allows much larger time evolution step size and reduces memory consumption. However, its setup is a little more complex than propagation and compression, since one-site TDVP generally can not adjust bond dimension during time evolution.
[22]:
from renormalizer.utils.configs import EvolveConfig, EvolveMethod
[23]:
mps.evolve_config = EvolveConfig(EvolveMethod.tdvp_ps)
Next, we perform one more step of the time evolution using TDVP-PS.
Note that the bond dimension in the middle of the MPS does not increase from 4 to 8. This is because one-site TDVP-PS does not alter the bond dimension.
In order to perform one-site TDVP-PS time evolution, it is imperative to increase the bond dimension of the MPS first, and then perform the time evolution.
[24]:
new_mps = mps.evolve(mpo, dt)
new_mps.bond_dims
[24]:
[1, 2, 4, 1]
The expand_bond_dimension
function increases the bond dimension to target value specified by the compress_config
. In short, the function adds the vectors in the Krylov space \(H^n|\psi\rangle\) to the input wavefunction \(|\psi\rangle\), and then compresses it to the target bond dimension.
The include_ex
option is specifically designed for systems with quantum number conservation. If the initial wavefunction does not include contributions from a certain symmetry sector, these contributions will not reappear during time evolution due to the projection error in TDVP. Setting include_ex=True
will add a small perturbation to the initial state to help recover the missing symmetry sector contributions. In our case, the quantum number conservation is disabled, so we set
include_ex=False
.
[25]:
new_mps = mps.expand_bond_dimension(mpo, include_ex=False)
new_mps.bond_dims
[25]:
[1, 2, 8, 1]
Observing Marcus Inverted Region#
Next, we put everything together and perform time evolution with different \(\Delta G\) using TDVP-PS. Since the initial state of time evolution is a Hartree product state, we must expand the bond dimension before performing the time evolution.
[26]:
# time evolution step
dt = 0.5
# initialize the MPS in the Hartree product state
init_mps = Mps.hartree_product_state(model, condition=init_condition)
logger.info(f"Initial MPS bond dimension: {init_mps.bond_dims}")
# setup compression configuration
init_mps.compress_config = CompressConfig(CompressCriteria.fixed, max_bonddim=8)
# setup time evolution configuration
init_mps.evolve_config = EvolveConfig(EvolveMethod.tdvp_ps)
# record the electronic occupation
n_list_list = []
delta_g_list = np.linspace(0, -2, 5)
for delta_g in delta_g_list:
# reconstruct the Hamiltonian. We can reuse the occupation operator though
model = get_model(delta_g)
mpo = Mpo(model)
# since using TDVP-PS, expand the bond dimension to target value.
# otherwise the time evolution is limited to the Hartree product state.
mps = init_mps.expand_bond_dimension(mpo, include_ex=False)
logger.info(f"MPS bond dimension after expanding: {mps.bond_dims}")
n_list = []
for i_step in range(20):
n = mps.expectation(n_op)
logger.info(f"Step {i_step}. Time {i_step * dt:.2f}. $n$ {n}")
# perform time evolution. Note that the evolution is not in-place.
mps = mps.evolve(mpo, dt)
n_list.append(n)
n_list_list.append(n_list)
2025-08-18 06:51:03,579[INFO] Initial MPS bond dimension: [1, 1, 1, 1]
2025-08-18 06:51:03,601[INFO] MPS bond dimension after expanding: [1, 2, 8, 1]
2025-08-18 06:51:03,602[INFO] Step 0. Time 0.00. $n$ 1.0
2025-08-18 06:51:03,629[INFO] Step 1. Time 0.50. $n$ 0.9975776788496079
2025-08-18 06:51:03,656[INFO] Step 2. Time 1.00. $n$ 0.9911341778211974
2025-08-18 06:51:03,682[INFO] Step 3. Time 1.50. $n$ 0.9824981838275648
2025-08-18 06:51:03,709[INFO] Step 4. Time 2.00. $n$ 0.9733283418283656
2025-08-18 06:51:03,736[INFO] Step 5. Time 2.50. $n$ 0.9645218516469746
2025-08-18 06:51:03,762[INFO] Step 6. Time 3.00. $n$ 0.9563382368492221
2025-08-18 06:51:03,789[INFO] Step 7. Time 3.50. $n$ 0.948730433389221
2025-08-18 06:51:03,816[INFO] Step 8. Time 4.00. $n$ 0.9415533981894397
2025-08-18 06:51:03,842[INFO] Step 9. Time 4.50. $n$ 0.9346400325046196
2025-08-18 06:51:03,869[INFO] Step 10. Time 5.00. $n$ 0.9278351842220715
2025-08-18 06:51:03,895[INFO] Step 11. Time 5.50. $n$ 0.921029899616409
2025-08-18 06:51:03,922[INFO] Step 12. Time 6.00. $n$ 0.9141749831330702
2025-08-18 06:51:03,948[INFO] Step 13. Time 6.50. $n$ 0.9072563582422936
2025-08-18 06:51:03,975[INFO] Step 14. Time 7.00. $n$ 0.900259128610462
2025-08-18 06:51:04,002[INFO] Step 15. Time 7.50. $n$ 0.8931521085342414
2025-08-18 06:51:04,030[INFO] Step 16. Time 8.00. $n$ 0.8858943878575144
2025-08-18 06:51:04,056[INFO] Step 17. Time 8.50. $n$ 0.8784553598367529
2025-08-18 06:51:04,084[INFO] Step 18. Time 9.00. $n$ 0.8708472420517187
2025-08-18 06:51:04,110[INFO] Step 19. Time 9.50. $n$ 0.8631616549410932
2025-08-18 06:51:04,157[INFO] MPS bond dimension after expanding: [1, 2, 8, 1]
2025-08-18 06:51:04,158[INFO] Step 0. Time 0.00. $n$ 1.0000000000000002
2025-08-18 06:51:04,185[INFO] Step 1. Time 0.50. $n$ 0.9975401925879046
2025-08-18 06:51:04,212[INFO] Step 2. Time 1.00. $n$ 0.9906013628490234
2025-08-18 06:51:04,238[INFO] Step 3. Time 1.50. $n$ 0.9802488858640296
2025-08-18 06:51:04,264[INFO] Step 4. Time 2.00. $n$ 0.9676593608104593
2025-08-18 06:51:04,290[INFO] Step 5. Time 2.50. $n$ 0.9537696503060078
2025-08-18 06:51:04,316[INFO] Step 6. Time 3.00. $n$ 0.9392072948520962
2025-08-18 06:51:04,343[INFO] Step 7. Time 3.50. $n$ 0.9243586104260211
2025-08-18 06:51:04,369[INFO] Step 8. Time 4.00. $n$ 0.9094397848359158
2025-08-18 06:51:04,396[INFO] Step 9. Time 4.50. $n$ 0.894551123228381
2025-08-18 06:51:04,422[INFO] Step 10. Time 5.00. $n$ 0.879743003866583
2025-08-18 06:51:04,448[INFO] Step 11. Time 5.50. $n$ 0.8650797798247613
2025-08-18 06:51:04,475[INFO] Step 12. Time 6.00. $n$ 0.8506519218071086
2025-08-18 06:51:04,501[INFO] Step 13. Time 6.50. $n$ 0.8365374567180583
2025-08-18 06:51:04,528[INFO] Step 14. Time 7.00. $n$ 0.8227688239625943
2025-08-18 06:51:04,554[INFO] Step 15. Time 7.50. $n$ 0.8093300617593504
2025-08-18 06:51:04,580[INFO] Step 16. Time 8.00. $n$ 0.796163296069385
2025-08-18 06:51:04,606[INFO] Step 17. Time 8.50. $n$ 0.7831685752085114
2025-08-18 06:51:04,632[INFO] Step 18. Time 9.00. $n$ 0.7702000677823149
2025-08-18 06:51:04,657[INFO] Step 19. Time 9.50. $n$ 0.7570646850987612
2025-08-18 06:51:04,703[INFO] MPS bond dimension after expanding: [1, 2, 8, 1]
2025-08-18 06:51:04,705[INFO] Step 0. Time 0.00. $n$ 1.0
2025-08-18 06:51:04,730[INFO] Step 1. Time 0.50. $n$ 0.9975277339592452
2025-08-18 06:51:04,756[INFO] Step 2. Time 1.00. $n$ 0.9904261839751384
2025-08-18 06:51:04,782[INFO] Step 3. Time 1.50. $n$ 0.9795268047068679
2025-08-18 06:51:04,808[INFO] Step 4. Time 2.00. $n$ 0.9659170449060821
2025-08-18 06:51:04,834[INFO] Step 5. Time 2.50. $n$ 0.9506850089406198
2025-08-18 06:51:04,859[INFO] Step 6. Time 3.00. $n$ 0.9347455146607804
2025-08-18 06:51:04,885[INFO] Step 7. Time 3.50. $n$ 0.9187422726641743
2025-08-18 06:51:04,910[INFO] Step 8. Time 4.00. $n$ 0.9030158921696072
2025-08-18 06:51:04,936[INFO] Step 9. Time 4.50. $n$ 0.8876520609099678
2025-08-18 06:51:04,962[INFO] Step 10. Time 5.00. $n$ 0.872598279802041
2025-08-18 06:51:04,988[INFO] Step 11. Time 5.50. $n$ 0.8577742342547119
2025-08-18 06:51:05,014[INFO] Step 12. Time 6.00. $n$ 0.8431151893156004
2025-08-18 06:51:05,040[INFO] Step 13. Time 6.50. $n$ 0.8285809886892017
2025-08-18 06:51:05,066[INFO] Step 14. Time 7.00. $n$ 0.8141772509631202
2025-08-18 06:51:05,092[INFO] Step 15. Time 7.50. $n$ 0.7999687732355262
2025-08-18 06:51:05,118[INFO] Step 16. Time 8.00. $n$ 0.7860556760297761
2025-08-18 06:51:05,142[INFO] Step 17. Time 8.50. $n$ 0.772520337668073
2025-08-18 06:51:05,168[INFO] Step 18. Time 9.00. $n$ 0.7593659426880239
2025-08-18 06:51:05,193[INFO] Step 19. Time 9.50. $n$ 0.7464615604976066
2025-08-18 06:51:05,238[INFO] MPS bond dimension after expanding: [1, 2, 8, 1]
2025-08-18 06:51:05,239[INFO] Step 0. Time 0.00. $n$ 1.0
2025-08-18 06:51:05,265[INFO] Step 1. Time 0.50. $n$ 0.9975406090105958
2025-08-18 06:51:05,290[INFO] Step 2. Time 1.00. $n$ 0.990624961140675
2025-08-18 06:51:05,316[INFO] Step 3. Time 1.50. $n$ 0.9804700824710465
2025-08-18 06:51:05,341[INFO] Step 4. Time 2.00. $n$ 0.9686190882357233
2025-08-18 06:51:05,367[INFO] Step 5. Time 2.50. $n$ 0.9564506630039298
2025-08-18 06:51:05,392[INFO] Step 6. Time 3.00. $n$ 0.9448244017405041
2025-08-18 06:51:05,418[INFO] Step 7. Time 3.50. $n$ 0.9339798618065525
2025-08-18 06:51:05,443[INFO] Step 8. Time 4.00. $n$ 0.9236970054402205
2025-08-18 06:51:05,468[INFO] Step 9. Time 4.50. $n$ 0.913611841996374
2025-08-18 06:51:05,493[INFO] Step 10. Time 5.00. $n$ 0.9034945341198809
2025-08-18 06:51:05,518[INFO] Step 11. Time 5.50. $n$ 0.8933290940558755
2025-08-18 06:51:05,543[INFO] Step 12. Time 6.00. $n$ 0.8832123669566792
2025-08-18 06:51:05,568[INFO] Step 13. Time 6.50. $n$ 0.8732206640791189
2025-08-18 06:51:05,593[INFO] Step 14. Time 7.00. $n$ 0.863341431359627
2025-08-18 06:51:05,619[INFO] Step 15. Time 7.50. $n$ 0.8534881635602912
2025-08-18 06:51:05,645[INFO] Step 16. Time 8.00. $n$ 0.8435889215797623
2025-08-18 06:51:05,670[INFO] Step 17. Time 8.50. $n$ 0.833692245809423
2025-08-18 06:51:05,694[INFO] Step 18. Time 9.00. $n$ 0.8239993933202343
2025-08-18 06:51:05,719[INFO] Step 19. Time 9.50. $n$ 0.8147658344266653
2025-08-18 06:51:05,764[INFO] MPS bond dimension after expanding: [1, 2, 8, 1]
2025-08-18 06:51:05,765[INFO] Step 0. Time 0.00. $n$ 1.0000000000000002
2025-08-18 06:51:05,790[INFO] Step 1. Time 0.50. $n$ 0.9975784978846655
2025-08-18 06:51:05,816[INFO] Step 2. Time 1.00. $n$ 0.9911783314788769
2025-08-18 06:51:05,842[INFO] Step 3. Time 1.50. $n$ 0.9828793267407632
2025-08-18 06:51:05,867[INFO] Step 4. Time 2.00. $n$ 0.974804737510615
2025-08-18 06:51:05,893[INFO] Step 5. Time 2.50. $n$ 0.9681008008243925
2025-08-18 06:51:05,919[INFO] Step 6. Time 3.00. $n$ 0.9626946118364343
2025-08-18 06:51:05,945[INFO] Step 7. Time 3.50. $n$ 0.9578334887382822
2025-08-18 06:51:05,971[INFO] Step 8. Time 4.00. $n$ 0.952863688361241
2025-08-18 06:51:05,996[INFO] Step 9. Time 4.50. $n$ 0.9476402488995228
2025-08-18 06:51:06,021[INFO] Step 10. Time 5.00. $n$ 0.942371478319725
2025-08-18 06:51:06,046[INFO] Step 11. Time 5.50. $n$ 0.9372226613912158
2025-08-18 06:51:06,071[INFO] Step 12. Time 6.00. $n$ 0.9321434873912731
2025-08-18 06:51:06,096[INFO] Step 13. Time 6.50. $n$ 0.9270168710222633
2025-08-18 06:51:06,121[INFO] Step 14. Time 7.00. $n$ 0.9218472771188858
2025-08-18 06:51:06,145[INFO] Step 15. Time 7.50. $n$ 0.9167421622570097
2025-08-18 06:51:06,170[INFO] Step 16. Time 8.00. $n$ 0.9117276598161915
2025-08-18 06:51:06,195[INFO] Step 17. Time 8.50. $n$ 0.9066479717963585
2025-08-18 06:51:06,221[INFO] Step 18. Time 9.00. $n$ 0.9013355144433687
2025-08-18 06:51:06,246[INFO] Step 19. Time 9.50. $n$ 0.8959335401185693
By plotting the figure, we see that when \(\Delta G=-1\), the charge transfer rate is highest.
[27]:
t = np.arange(len(n_list_list[0])) * dt
for i, n_list in enumerate(n_list_list):
label = r"$\Delta G =" + str(delta_g_list[i]) + "$"
plt.plot(t, n_list, label=label, linestyle="--")
plt.xlabel("$t$")
plt.ylabel(r"Occupation $\langle a^\dagger_0 a_0 \rangle$")
plt.legend()
[27]:
<matplotlib.legend.Legend at 0x7fb327bcfc10>

The prediction is consistent with the Marcus theory ,where the Marcus rate is
Here \(\lambda = 2g^2\omega = 1\) is the reorganization energy.
Features of Time Evolution Schemes#
Renormalizer has a lot of different time evolution schemes, the most commonly used ones are 1. Propagation and Compression: EvolveMethod.prop_and_compress
2. Single-site TDVP-PS: EvolveMethod.tdvp_ps
3. Two-site TDVP-PS: EvolveMethod.tdvp_ps2
4. TDVP with varianble mean field (VMF): EvolveMethod.tdvp_vmf
The features of these methods are as follows
Name |
Accuracy |
Computational Cost |
Adaptive Step-Size |
Adaptive Bond Dimension |
Projection Error |
---|---|---|---|---|---|
|
Low |
High |
partial |
Yes |
No |
|
High |
Low |
No |
No |
Yes |
|
High |
High |
No |
Yes |
Partial |
|
High |
High |
Yes |
No |
Yes |
“Adaptive Step-Size” refers to the ability to dynamically adjust the time-evolution step size based on a predefined error threshold. In practice, this means the method automatically takes smaller steps when the wavefunction changes abruptly and larger steps when it evolves smoothly. For example, tdvp_vmf
performs multiple small steps within a single call to the evolve
method. This capability removes the need for users to manually check step-size convergence in TD-DMRG by testing multiple
step sizes.
“Adaptive Bond Dimension” refers to the ability to dynamically adjust the bond dimension during time evolution. Typically, the initial state is a Hartree product state. For tdvp_ps
and tdvp_vmf
, the expand_bond_dimension
method should be called before starting the time evolution. Otherwise, the simulation will proceed with a bond dimension of 1 rather than the value specified in the optimization settings.
“Projection Error” is the error introduced when projecting the time derivative onto the MPS manifold.
Empirically, we recommend tdvp_ps
for production level calculations.
[ ]: