Low Level Data Structure and Algorithms#
Matrix Product (MP)#
- class renormalizer.mps.mp.MatrixProduct[source]#
-
- _get_big_qn(cidx: List[int], swap=False)[source]#
get the quantum number of L-block and R-block renormalized basis
- Parameters:
cidx (list) – a list of center(active) site index. For 1site/2site algorithm, cidx has one/two elements.
- Returns:
qnbigl (np.ndarray) – super-L-block (L-block + active site if necessary) quantum number
qnbigr (np.ndarray) – super-R-block (active site + R-block if necessary) quantum number
qnmat (np.ndarray) – L-block + active site + R-block quantum number
- _update_mps(cstruct, cidx, qnbigl, qnbigr, percent=0)[source]#
update mps with basis selection algorithm of J. Chem. Phys. 120, 3172 (2004).
- Parameters:
cstruct (ndarray, List[ndarray]) – The active site coefficient.
cidx (list) – The List of active site index.
qnbigl (ndarray) – The super-L-block quantum number.
qnbigr (ndarray) – The super-R-block quantum number.
percent (float, int) – The percentage of renormalized basis which is equally selected from each quantum number section rather than according to singular values.
percent
is defined inprocedure
ofrenormalizer.utils.configs.OptimizeConfig
andvprocedure
ofrenormalizer.utils.configs.CompressConfig
.
- Returns:
- if
cstruct
is a list,averaged_ms
is a list of rotated ms of each element in
cstruct
as a single site calculation. It is used for better initial guess in SA-DMRG algorithm. Otherwise,None
is returned.self
is overwritten inplace.
- if
- Return type:
averaged_ms
- _update_ms(idx: int, u: ndarray, vt: ndarray, sigma=None, qnlset=None, qnrset=None, m_trunc=None)[source]#
update mps directly after svd
- add(other: MatrixProduct)[source]#
- check_left_canonical(rtol: Optional[float] = None, atol: Optional[float] = None)[source]#
check L-canonical
- check_right_canonical(rtol: Optional[float] = None, atol: Optional[float] = None)[source]#
check R-canonical
- compress(temp_m_trunc=None, ret_s=False)[source]#
inp: canonicalise MPS (or MPO)
- Return type:
truncated MPS
- dot(other: MatrixProduct) complex [source]#
dot product of two mps / mpo
- dot_ob(other: MatrixProduct) complex [source]#
dot product of two mps / mpo with open boundary, but the boundary of mps/mpo is larger than 1, different from the normal mps/mpo
- property is_complex#
- property is_left_canonical#
check the qn center in the L-canonical structure
- property is_mpdm#
- property is_mpo#
- property is_mps#
- property is_right_canonical#
check the qn center in the R-canonical structure
- metacopy() MatrixProduct [source]#
- move_qnidx(dstidx: int)[source]#
Quantum number has a boundary site, left hand of the site is L system qn, right hand of the side is R system qn, the sum of quantum number of L system and R system is tot.
- property mp_norm: float#
- if self.is_left_canon:
assert self.check_left_canonical() return np.linalg.norm(np.ravel(self[-1]))
- else:
assert self.check_right_canonical() return np.linalg.norm(np.ravel(self[0]))
- property pbond_dims#
- property pbond_list#
- property site_num#
- property threshold#
- property total_bytes#
- variational_compress(mpo=None, guess=None)[source]#
Variational compress an mps/mpdm/mpo
- Parameters:
mpo (renormalizer.mps.Mpo, optional) – Default is
None
. if mpo is notNone
, the returned mps is an approximation ofmpo @ self
guess (renormalizer.mps.MatrixProduct, optional) – Initial guess of compressed mps/mpdm/mpo. Default is
None
.
Note
the variational compress related configurations is defined in
self
ifguess=None
, otherwise is defined inguess
- Returns:
mp – a new compressed mps/mpdm/mpo,
self
is not overwritten.guess
is overwritten.- Return type:
renormalizer.mps.MatrixProduct
Matrix Product State (MPS)#
- class renormalizer.mps.Mps[source]#
-
- angle(other)#
- append(array)#
- build_empty_mp(num)#
- build_empty_qn()#
- build_none_qn()#
- calc_1site_rdm(idx=None)[source]#
Calculate 1-site reduced density matrix
\(\rho_i = \textrm{Tr}_{j \neq i} | \Psi \rangle \langle \Psi|\)
- calc_2site_mutual_entropy() ndarray [source]#
Calculate mutual entropy between two sites. Also known as mutual information
\(m_{ij} = (s_i + s_j - s_{ij})/2\)
See Chemical Physics 323 (2006) 519–531
- Returns:
mutual_entropy – mutual entropy with shape (nsite, nsite)
- Return type:
2d np.ndarry
- calc_2site_rdm()[source]#
Calculate 2-site reduced density matrix
\(\rho_{ij} = \textrm{Tr}_{k \neq i, k \neq j} | \Psi \rangle \langle \Psi |\).
- Returns:
rdm – \(\{(0,1):\rho_{01}, (0,2):\rho_{02}, \cdots\}\). The key is a tuple of index of the site.
- Return type:
Dict
- calc_bond_entropy() ndarray [source]#
Calculate von Neumann entropy at each bond according to \(S = -\textrm{Tr}(\rho \ln \rho)\) where \(\rho\) is the reduced density matrix of either block.
- Returns:
S – a NumPy array containing the entropy values.
- Return type:
1D array
- calc_edof_rdm() ndarray [source]#
Calculate the reduced density matrix of electronic DoF
\(\rho_{ij} = \langle \Psi | a_i^\dagger a_j | \Psi \rangle\)
Note
This function is designed for single-electron system. Fermionic commutation relation is not considered.
- calc_entropy(entropy_type)[source]#
Calculate 1site, 2site, mutual and bond Von Neumann entropy
\(\textrm{entropy} = -\textrm{Tr}(\rho \ln \rho)\) where \(\ln\) stands for natural logarithm.
1site entropy is the entropy between any site and the other
(N-1)
sites. 2site entropy is the entropy between any two sites and the other(N-2)
sites. mutual entropy characterize the entropy between any two sites. bond entropy is the entropy between L-block and R-block.
- check_right_canonical(rtol: Optional[float] = None, atol: Optional[float] = None)#
check R-canonical
- compress(temp_m_trunc=None, ret_s=False)#
inp: canonicalise MPS (or MPO)
- Return type:
truncated MPS
- copy()#
- property digest#
- dot(other: MatrixProduct) complex #
dot product of two mps / mpo
- dot_ob(other: MatrixProduct) complex #
dot product of two mps / mpo with open boundary, but the boundary of mps/mpo is larger than 1, different from the normal mps/mpo
- property e_occupations#
Electronic occupations \(a^\dagger_i a_i\) for each electronic DoF. The order is defined by
e_dofs
.
- expand_bond_dimension(hint_mpo=None, coef=1e-10, include_ex=True)[source]#
expand bond dimension as required in compress_config
- classmethod from_mp(model, mplist)#
- classmethod ground_state(model: Model, max_entangled: bool, normalize: bool = True, condition: Optional[Dict] = None)[source]#
Obtain ground state at \(T = 0\) or \(T = \infty\) (maximum entangled). Electronic DOFs are always at ground state. and vibrational DOFs depend on
max_entangled
. For Spin-Boson model the electronic DOF also depends onmax_entangled
.- Parameters:
model (
Model
) – system information.max_entangled (bool) – temperature of the vibrational DOFs. If set to
True
, \(T = \infty\) and if set toFalse
, \(T = 0\).normalize (bool, optional) – if the returned Mps are normalized when
max_entangled=True
. Default is True. Ifnormalize=False
, the vibrational part is identity.condition (dict, optional) – the same as
hartree_product_state
. only used in ba.BasisMultiElectron cases to define the occupation. Default isNone
.
- Returns:
mps
- Return type:
- classmethod hartree_product_state(model, condition: Dict, qn_idx: Optional[int] = None)[source]#
Construct a Hartree product state
- Parameters:
model (
Model
) – Model information.condition (Dict) –
Dict with format
{dof:local_state}
. The default local state for dofs not specified is the “0” state. An example is{"e_1":1, "v_0":2, "v_3":[0, 0.707, 0.707]}
.Note
If there are bases that contain multiple dofs in the model, the value of the dict is the state of all dofs of the basis. For example, if a basis contains
"e_1"
,"e_2"
and"e_3"
,{"e_1": 2}
({"e_1": [0, 0, 1]}
) means"e_3"
is occupied and{"e_1": 1}
({"e_1": [0, 1, 0]}
) means"e_2"
is occupied. Be aware that inrenormalizer.BasisMultiElectronVac
the vacuum state is added to the0
index.qn_idx (int) – the site index of the quantum number center.
- Returns:
Constructed mps (
Mps
)
- property is_complex#
- property is_left_canonical#
check the qn center in the L-canonical structure
- property is_mpdm#
- property is_mpo#
- property is_mps#
- property is_right_canonical#
check the qn center in the R-canonical structure
- move_qnidx(dstidx: int)#
Quantum number has a boundary site, left hand of the site is L system qn, right hand of the side is R system qn, the sum of quantum number of L system and R system is tot.
- property mp_norm: float#
- if self.is_left_canon:
assert self.check_left_canonical() return np.linalg.norm(np.ravel(self[-1]))
- else:
assert self.check_right_canonical() return np.linalg.norm(np.ravel(self[0]))
- property nexciton#
- property norm#
the norm of the total wavefunction
- normalize(kind)[source]#
normalize the wavefunction
- Parameters:
kind (str) – “mps_only”: the mps part is normalized and coeff is not modified; “mps_norm_to_coeff”: the mps part is normalized and the norm is multiplied to coeff; “mps_and_coeff”: both mps and coeff is normalized
- Return type:
self
is overwritten.
- property pbond_dims#
- property pbond_list#
- property ph_occupations#
phonon occupations \(b^\dagger_i b_i\) for each electronic DoF. The order is defined by
v_dofs
.
- scale(val, inplace=False)#
- property site_num#
- property threshold#
- property total_bytes#
- variational_compress(mpo=None, guess=None)#
Variational compress an mps/mpdm/mpo
- Parameters:
mpo (renormalizer.mps.Mpo, optional) – Default is
None
. if mpo is notNone
, the returned mps is an approximation ofmpo @ self
guess (renormalizer.mps.MatrixProduct, optional) – Initial guess of compressed mps/mpdm/mpo. Default is
None
.
Note
the variational compress related configurations is defined in
self
ifguess=None
, otherwise is defined inguess
- Returns:
mp – a new compressed mps/mpdm/mpo,
self
is not overwritten.guess
is overwritten.- Return type:
renormalizer.mps.MatrixProduct
Matrix Product Operator (MPO)#
- class renormalizer.mps.Mpo(model: ~typing.Optional[~renormalizer.model.model.Model] = None, terms: ~typing.Optional[~typing.Union[~renormalizer.model.op.Op, ~typing.List[~renormalizer.model.op.Op]]] = None, offset: ~renormalizer.utils.quantity.Quantity = <renormalizer.utils.quantity.Quantity object>, algo='qr')[source]#
Matrix product operator (MPO)
- add(other: MatrixProduct)#
- angle(other)#
- append(array)#
- apply(mp: MatrixProduct, canonicalise: bool = False) MatrixProduct [source]#
- build_empty_mp(num)#
- build_empty_qn()#
- build_none_qn()#
- check_right_canonical(rtol: Optional[float] = None, atol: Optional[float] = None)#
check R-canonical
- compress(temp_m_trunc=None, ret_s=False)#
inp: canonicalise MPS (or MPO)
- Return type:
truncated MPS
- conj()#
complex conjugate
- contract(mps, algo='svd')[source]#
an approximation of mpo @ mps/mpdm/mpo
- Parameters:
- Returns:
new_mps – an approximation of mpo @ mps/mpdm/mpo. The input
mps
is not overwritten.- Return type:
See also
renormalizer.mps.mp.MatrixProduct.compress
svd compression.
renormalizer.mps.mp.MatrixProduct.variational_compress
variational compression.
- copy()#
- property digest#
- dot(other: MatrixProduct) complex #
dot product of two mps / mpo
- dot_ob(other: MatrixProduct) complex #
dot product of two mps / mpo with open boundary, but the boundary of mps/mpo is larger than 1, different from the normal mps/mpo
- property dummy_qn#
- dump(fname, other_attrs=None)#
- classmethod exact_propagator(model: HolsteinModel, x, space='GS', shift=0.0)[source]#
construct the GS space propagator e^{xH} exact MPO H=sum_{in} omega_{in} b^dagger_{in} b_{in} fortunately, the H is local. so e^{xH} = e^{xh1}e^{xh2}…e^{xhn} the bond dimension is 1 shift is the a constant for H+shift
- classmethod from_mp(model, mplist)#
- classmethod intersite(model: ~renormalizer.model.model.HolsteinModel, e_opera: dict, ph_opera: dict, scale: ~renormalizer.utils.quantity.Quantity = <renormalizer.utils.quantity.Quantity object>)[source]#
construct the inter site MPO
- Parameters:
model (HolsteinModel) – the molecular information
e_opera – the electronic operators. {imol: operator}, such as {1:”a”, 3:r”a^dagger”}
ph_opera – the vibrational operators. {(imol, iph): operator}, such as {(0,5):”b”}
scale (Quantity) – scalar to scale the mpo
Note
the operator index starts from 0,1,2…
- property is_complex#
- property is_left_canonical#
check the qn center in the L-canonical structure
- property is_mpdm#
- property is_mpo#
- property is_mps#
- property is_right_canonical#
check the qn center in the R-canonical structure
- move_qnidx(dstidx: int)#
Quantum number has a boundary site, left hand of the site is L system qn, right hand of the side is R system qn, the sum of quantum number of L system and R system is tot.
- property mp_norm: float#
- if self.is_left_canon:
assert self.check_left_canonical() return np.linalg.norm(np.ravel(self[-1]))
- else:
assert self.check_right_canonical() return np.linalg.norm(np.ravel(self[0]))
- property pbond_dims#
- property pbond_list#
- classmethod ph_onsite(model: HolsteinModel, opera: str, mol_idx: int, ph_idx=0)[source]#
- scale(val, inplace=False)#
- property site_num#
- property threshold#
- to_complex(inplace=False)#
- property total_bytes#
- variational_compress(mpo=None, guess=None)#
Variational compress an mps/mpdm/mpo
- Parameters:
mpo (renormalizer.mps.Mpo, optional) – Default is
None
. if mpo is notNone
, the returned mps is an approximation ofmpo @ self
guess (renormalizer.mps.MatrixProduct, optional) – Initial guess of compressed mps/mpdm/mpo. Default is
None
.
Note
the variational compress related configurations is defined in
self
ifguess=None
, otherwise is defined inguess
- Returns:
mp – a new compressed mps/mpdm/mpo,
self
is not overwritten.guess
is overwritten.- Return type:
renormalizer.mps.MatrixProduct
Matrix Product Density Matrix (MPDM)#
- class renormalizer.mps.MpDm[source]#
- add(other)#
- angle(other)#
- append(array)#
- build_empty_mp(num)#
- build_empty_qn()#
- build_none_qn()#
- calc_1site_rdm(idx=None)#
Calculate 1-site reduced density matrix
\(\rho_i = \textrm{Tr}_{j \neq i} | \Psi \rangle \langle \Psi|\)
- calc_2site_mutual_entropy() ndarray #
Calculate mutual entropy between two sites. Also known as mutual information
\(m_{ij} = (s_i + s_j - s_{ij})/2\)
See Chemical Physics 323 (2006) 519–531
- Returns:
mutual_entropy – mutual entropy with shape (nsite, nsite)
- Return type:
2d np.ndarry
- calc_2site_rdm()#
Calculate 2-site reduced density matrix
\(\rho_{ij} = \textrm{Tr}_{k \neq i, k \neq j} | \Psi \rangle \langle \Psi |\).
- Returns:
rdm – \(\{(0,1):\rho_{01}, (0,2):\rho_{02}, \cdots\}\). The key is a tuple of index of the site.
- Return type:
Dict
- calc_bond_entropy() ndarray #
Calculate von Neumann entropy at each bond according to \(S = -\textrm{Tr}(\rho \ln \rho)\) where \(\rho\) is the reduced density matrix of either block.
- Returns:
S – a NumPy array containing the entropy values.
- Return type:
1D array
- calc_edof_rdm() ndarray #
Calculate the reduced density matrix of electronic DoF
\(\rho_{ij} = \langle \Psi | a_i^\dagger a_j | \Psi \rangle\)
Note
This function is designed for single-electron system. Fermionic commutation relation is not considered.
- calc_entropy(entropy_type)#
Calculate 1site, 2site, mutual and bond Von Neumann entropy
\(\textrm{entropy} = -\textrm{Tr}(\rho \ln \rho)\) where \(\ln\) stands for natural logarithm.
1site entropy is the entropy between any site and the other
(N-1)
sites. 2site entropy is the entropy between any two sites and the other(N-2)
sites. mutual entropy characterize the entropy between any two sites. bond entropy is the entropy between L-block and R-block.
- check_right_canonical(rtol: Optional[float] = None, atol: Optional[float] = None)#
check R-canonical
- compress(temp_m_trunc=None, ret_s=False)#
inp: canonicalise MPS (or MPO)
- Return type:
truncated MPS
- contract(mps, algo='svd')#
an approximation of mpo @ mps/mpdm/mpo
- Parameters:
- Returns:
new_mps – an approximation of mpo @ mps/mpdm/mpo. The input
mps
is not overwritten.- Return type:
See also
renormalizer.mps.mp.MatrixProduct.compress
svd compression.
renormalizer.mps.mp.MatrixProduct.variational_compress
variational compression.
- copy()#
- property digest#
- dot(other: MatrixProduct) complex #
dot product of two mps / mpo
- dot_ob(other: MatrixProduct) complex #
dot product of two mps / mpo with open boundary, but the boundary of mps/mpo is larger than 1, different from the normal mps/mpo
- property dummy_qn#
- dump(fname)#
- property e_occupations#
Electronic occupations \(a^\dagger_i a_i\) for each electronic DoF. The order is defined by
e_dofs
.
- classmethod exact_propagator(model: HolsteinModel, x, space='GS', shift=0.0)#
construct the GS space propagator e^{xH} exact MPO H=sum_{in} omega_{in} b^dagger_{in} b_{in} fortunately, the H is local. so e^{xH} = e^{xh1}e^{xh2}…e^{xhn} the bond dimension is 1 shift is the a constant for H+shift
- expand_bond_dimension(hint_mpo=None, coef=1e-10, include_ex=True)#
expand bond dimension as required in compress_config
- expectations(mpos, self_conj=None, opt=True) ndarray #
- classmethod finiteT_cv(model, nexciton, m_max, spectratype, percent=1.0)#
- classmethod from_mp(model, mplist)#
- classmethod ground_state(model, max_entangled)[source]#
Obtain ground state at \(T = 0\) or \(T = \infty\) (maximum entangled). Electronic DOFs are always at ground state. and vibrational DOFs depend on
max_entangled
. For Spin-Boson model the electronic DOF also depends onmax_entangled
.- Parameters:
model (
Model
) – system information.max_entangled (bool) – temperature of the vibrational DOFs. If set to
True
, \(T = \infty\) and if set toFalse
, \(T = 0\).normalize (bool, optional) – if the returned Mps are normalized when
max_entangled=True
. Default is True. Ifnormalize=False
, the vibrational part is identity.condition (dict, optional) – the same as
hartree_product_state
. only used in ba.BasisMultiElectron cases to define the occupation. Default isNone
.
- Returns:
mps
- Return type:
- classmethod hartree_product_state(model, condition: Dict, qn_idx: Optional[int] = None)#
Construct a Hartree product state
- Parameters:
model (
Model
) – Model information.condition (Dict) –
Dict with format
{dof:local_state}
. The default local state for dofs not specified is the “0” state. An example is{"e_1":1, "v_0":2, "v_3":[0, 0.707, 0.707]}
.Note
If there are bases that contain multiple dofs in the model, the value of the dict is the state of all dofs of the basis. For example, if a basis contains
"e_1"
,"e_2"
and"e_3"
,{"e_1": 2}
({"e_1": [0, 0, 1]}
) means"e_3"
is occupied and{"e_1": 1}
({"e_1": [0, 1, 0]}
) means"e_2"
is occupied. Be aware that inrenormalizer.BasisMultiElectronVac
the vacuum state is added to the0
index.qn_idx (int) – the site index of the quantum number center.
- Returns:
Constructed mps (
Mps
)
- classmethod intersite(model: ~renormalizer.model.model.HolsteinModel, e_opera: dict, ph_opera: dict, scale: ~renormalizer.utils.quantity.Quantity = <renormalizer.utils.quantity.Quantity object>)#
construct the inter site MPO
- Parameters:
model (HolsteinModel) – the molecular information
e_opera – the electronic operators. {imol: operator}, such as {1:”a”, 3:r”a^dagger”}
ph_opera – the vibrational operators. {(imol, iph): operator}, such as {(0,5):”b”}
scale (Quantity) – scalar to scale the mpo
Note
the operator index starts from 0,1,2…
- property is_complex#
- is_hermitian()#
- property is_left_canonical#
check the qn center in the L-canonical structure
- property is_mpdm#
- property is_mpo#
- property is_mps#
- property is_right_canonical#
check the qn center in the R-canonical structure
- classmethod max_entangled_ex(model, normalize=True)[source]#
T = infty locally maximal entangled EX state
- move_qnidx(dstidx: int)#
Quantum number has a boundary site, left hand of the site is L system qn, right hand of the side is R system qn, the sum of quantum number of L system and R system is tot.
- property mp_norm: float#
- if self.is_left_canon:
assert self.check_left_canonical() return np.linalg.norm(np.ravel(self[-1]))
- else:
assert self.check_right_canonical() return np.linalg.norm(np.ravel(self[0]))
- property nexciton#
- property norm#
the norm of the total wavefunction
- normalize(kind)#
normalize the wavefunction
- Parameters:
kind (str) – “mps_only”: the mps part is normalized and coeff is not modified; “mps_norm_to_coeff”: the mps part is normalized and the norm is multiplied to coeff; “mps_and_coeff”: both mps and coeff is normalized
- Return type:
self
is overwritten.
- property pbond_dims#
- property pbond_list#
- property ph_occupations#
phonon occupations \(b^\dagger_i b_i\) for each electronic DoF. The order is defined by
v_dofs
.
- classmethod ph_onsite(model: HolsteinModel, opera: str, mol_idx: int, ph_idx=0)#
- promote_mt_type(mp)#
- scale(val, inplace=False)#
- property site_num#
- property threshold#
- property total_bytes#
- variational_compress(mpo=None, guess=None)#
Variational compress an mps/mpdm/mpo
- Parameters:
mpo (renormalizer.mps.Mpo, optional) – Default is
None
. if mpo is notNone
, the returned mps is an approximation ofmpo @ self
guess (renormalizer.mps.MatrixProduct, optional) – Initial guess of compressed mps/mpdm/mpo. Default is
None
.
Note
the variational compress related configurations is defined in
self
ifguess=None
, otherwise is defined inguess
- Returns:
mp – a new compressed mps/mpdm/mpo,
self
is not overwritten.guess
is overwritten.- Return type:
renormalizer.mps.MatrixProduct
Thermal Propagation#
- class renormalizer.mps.thermalprop.ThermalProp(init_mpdm: MpDm, h_mpo_model: Optional[Model] = None, exact: bool = False, space: str = 'GS', evolve_config: Optional[EvolveConfig] = None, dump_mps: Optional[bool] = None, dump_dir: Optional[str] = None, job_name: Optional[str] = None, properties: Optional[Property] = None, auto_expand: bool = True)[source]#
Bases:
TdMpsJob
Thermally Propagate initial matrix product density operator (
MpDm
) in imaginary time.- Parameters:
init_mpdm (
MpDm
) – the initial density matrix to be propagated. Usually identity.h_mpo_model (
Model
) – Model for the system Hamiltonian. Default is the same withinit_mpdm.model
.exact (bool) – whether propagate assuming Hamiltonian is local \(\hat H = \sum_i \hat H_i = \sum_{in} \omega_{in} b^\dagger_{in} b_{in}\) and exact propagation is possible through \(e^{xH} = e^{xh_1}e^{xh_2} \cdots e^{xh_n}\). If set to
True
, properties such as occupations are not calculated during time evolution for better efficiency.space (str) – the space of exact propagation. Possible options are
"GS"
or"EX"
. If set to"GS"
, then the exact propagation is performed in zero exciton space. If set to"EX"
, then the exact propagation is performed in one exciton space, i.e. the vibrations are regarded as displaced oscillators.evolve_config (
EvolveConfig
) – config when evolving the MpDm in imaginary time.dump_mps (str) – if dump mps when dumping, “all”, “one”, None; Default to None
dump_dir (str) – the directory for logging and numerical result output.
job_name (str) – the name of the calculation job which determines the file name of the logging and numerical result output.
properties (
Property
) –
- dump_dict()#
- property e_occupations_array#
- evolve(evolve_dt=None, nsteps=None, evolve_time=None)[source]#
- Parameters:
evolve_dt (float) – the time step to run
evolve_single_step
andprocess_mps
.nsteps (int) – the total number of evolution steps
evolve_time (float) – the total evolution time
Note
evolve_dt
math: ` imes`nsteps
=evolve_time
, otherwise nsteps has a higher priority.
- evolve_single_step(evolve_dt)[source]#
Evolve the mps for a single step with step size
evolve_dt
.- Returns:
new mps after the evolution
- property evolve_times_array#
- get_dump_dict()[source]#
Obtain calculated properties to dump in
dict
type.- Returns:
return a (ordered) dict to dump as npz
- property latest_evolve_time#
- property ph_occupations_array#
- process_mps(mps)[source]#
Process the newly evolved the mps. Primarily for the calculation of properties. Note that currently
self.latest_mps
has not been updated.- Parameters:
mps – The evolved new mps of the time step.
- stop_evolve_criteria()#
- property vn_entropy_array#
- renormalizer.mps.thermalprop.load_thermal_state(model, path: str)[source]#
Load thermal propagated state from disk. Return None if the file is not found.
- Parameters:
model (
MolList
) – system informationpath (str) – the path to load thermal state from. Should be an numpy
.npz
file.
Returns: Loaded MpDm