![]() |
Octopus
|
Data Types | |
type | dmp_t |
Density matrix dissipation. More... | |
Functions/Subroutines | |
subroutine, public | dm_propagation_init_run (adiabatic_st, namespace, space, gr, st, hm, mc) |
Initialise the adiabatic states prior to running TB propagation. More... | |
subroutine | dmp_init (this, namespace, st) |
Initialise an instance of density matrix dissipation. More... | |
subroutine, public | dm_propagation_run (dmp, namespace, space, gr, ions, st, mc, hm, ks, iter, dt, ext_partners, update_energy) |
Apply dissipation to a TD run via the Linblad formalism. More... | |
subroutine | orthogonality_check_ks (ik, st, gr) |
Check orthonality and electron number of TDKS wavefunctions We constructed the density matrix from the TDKS wavefunctions, as. More... | |
subroutine | construct_residuals (adiabatic_st, ik, st, gr, conj_psi_phi, resd) |
Construct the residual basis from the TDKS wavefunctions. More... | |
subroutine | construct_density_matrix (ik, st, conj_psi_phi, rho_mat) |
Construct the full density matrix in the adiabatic and residual basis. More... | |
subroutine | population_in_adiabatic (ik, st, rho_mat, pop) |
Calculate population in adiabatic basis. More... | |
subroutine | update_st (dmp, ik, gr, resd, st, rho_mat) |
Diagonalize the density matrix and update the wavefunction and occupation. More... | |
subroutine | update_wfc_occ (adiabatic_st, ik, st, gr, resd, rho_mat, occ) |
Transform the wavefunctions into real-space basis The wavefunctions in real-space basis are given by: More... | |
subroutine | update_wfc_occ_procrustes (adiabatic_st, ik, st, gr, resd, rho_mat, occ) |
Make TDKS wavefunctions continuous by maximizing their overlap. More... | |
subroutine | dissipation (dmp, bst, st, ik, dt, nn) |
Apply dissipation to the density matrix This subroutine applies the dissipation operator to the density matrix using Taylor expansion of the matrix exponential upto order 4. More... | |
subroutine, public dm_propagation_oct_m::dm_propagation_init_run | ( | type(states_elec_t), intent(inout) | adiabatic_st, |
type(namespace_t), intent(in) | namespace, | ||
type(electron_space_t), intent(in) | space, | ||
type(grid_t), intent(in) | gr, | ||
type(states_elec_t), intent(in) | st, | ||
type(hamiltonian_elec_t), intent(in) | hm, | ||
type(multicomm_t), intent(in) | mc | ||
) |
Initialise the adiabatic states prior to running TB propagation.
Definition at line 190 of file dm_propagation.F90.
|
private |
Initialise an instance of density matrix dissipation.
Definition at line 220 of file dm_propagation.F90.
subroutine, public dm_propagation_oct_m::dm_propagation_run | ( | type(dmp_t), intent(inout) | dmp, |
type(namespace_t), intent(in) | namespace, | ||
type(electron_space_t), intent(in) | space, | ||
type(grid_t), intent(in) | gr, | ||
type(ions_t), intent(in) | ions, | ||
type(states_elec_t), intent(inout) | st, | ||
type(multicomm_t), intent(in) | mc, | ||
type(hamiltonian_elec_t), intent(inout) | hm, | ||
type(v_ks_t), intent(inout) | ks, | ||
integer, intent(in) | iter, | ||
real(real64), intent(in) | dt, | ||
type(partner_list_t), intent(in) | ext_partners, | ||
logical, intent(in), optional | update_energy | ||
) |
Apply dissipation to a TD run via the Linblad formalism.
[in,out] | dmp | Density matrix propagation instance |
[in] | space | Spin dimension |
[in] | gr | Underlying mesh |
[in] | ions | Geometry |
[in,out] | st | Time-dependent wfcs |
[in] | mc | MPI communicator |
[in,out] | hm | Hamiltonian, updated at each time step |
[in,out] | ks | Kohn-Sham potential |
[in] | iter | Number of iterations |
[in] | dt | Time step |
[in] | ext_partners | External interaction partners |
[in] | update_energy | Whether update energy when solve adabatic wfcs |
Definition at line 313 of file dm_propagation.F90.
|
private |
Check orthonality and electron number of TDKS wavefunctions We constructed the density matrix from the TDKS wavefunctions, as.
\[ \rho_{ij} = \langle \psi_i | \left( \sum_l f_l | \psi_l \rangle \langle \psi_l | \right) | \psi_j \rangle \]
where \( f_l \) is the occupation of the state \(l\) at \(k\)-point.
[in] | ik | K-point index |
[in] | st | Time-dependent wfcs |
[in] | gr | Underlying mesh |
Definition at line 412 of file dm_propagation.F90.
|
private |
Construct the residual basis from the TDKS wavefunctions.
The residual basis vectors are defined as
\[ | d_n \rangle = | \psi_n \rangle - \sum_{i=1}^{N} | \phi_i \rangle \langle \phi_i | \psi_n \rangle \]
where \( | \psi_n \rangle \) are the TDKS wavefunctions, and \( | \phi_i \rangle \) are the adiabatic states.
The resulting residuals \( | d_n \rangle \) are stored in resd
, a 3D array with shape (np, dim, nst)
.
[in] | adiabatic_st | Adiabatic basis |
[in] | ik | K-point index |
[in] | st | Time-dependent wfcs |
[in] | gr | Underlying mesh |
[out] | conj_psi_phi | (i, j) = <\phi, d_j | \psi_i > |
[out] | resd | Residuals |
Definition at line 481 of file dm_propagation.F90.
|
private |
Construct the full density matrix in the adiabatic and residual basis.
The density matrix elements in the combined basis are given by
\[ \rho_{ij} = \langle \chi_i | \left( \sum_l f_l\, | \psi_l \rangle \langle \psi_l | \right) | \chi_j \rangle \]
where \( \{ \chi_i \} = \{ \phi_1, \dots, \phi_N, d_1, \dots, d_N \} \) is the complete basis, consisting of adiabatic states \( \phi_i \) and residuals \( d_i \), and \( f_l \) is the occupation of state \( \psi_l \).
[in] | ik | k-point index |
[in] | st | Time-dependent wfcs |
[in] | conj_psi_phi | (i, j) = <\phi, d_j | \psi_i > |
[out] | rho_mat | The full density matrix |
Definition at line 543 of file dm_propagation.F90.
|
private |
Calculate population in adiabatic basis.
Number of electrons in the adiabatic basis is
\[ pop = \sum_{i,k} \rho_{ii} w_k \]
[in] | ik | K-point index |
[in] | st | Time-dependent wfcs |
[in] | rho_mat | The full density matrix |
[in,out] | pop | Population in the adiabatic basis |
Definition at line 577 of file dm_propagation.F90.
|
private |
Diagonalize the density matrix and update the wavefunction and occupation.
This routine solves the generalized eigenvalue problem:
\[ \rho \psi = \lambda B \psi, \]
where \( \rho \) is the density matrix, \( B \) is the overlap matrix, and \( \psi \) are the wavefunctions.
Definition at line 608 of file dm_propagation.F90.
|
private |
Transform the wavefunctions into real-space basis The wavefunctions in real-space basis are given by:
\[ \psi = \sum_{i=1}^{N} c_i \phi_i + c_{i+N} \, \mathrm{resd}_i, \]
where \( \phi_i \) are the adiabatic states and \( \mathrm{resd}_i \) are the residual basis vectors.
Definition at line 685 of file dm_propagation.F90.
|
private |
Make TDKS wavefunctions continuous by maximizing their overlap.
Applies a unitary transformation to \( \tilde{\Psi} \) such that the overlap with the input wavefunctions is maximized.
The resulting wavefunctions \( \psi_i \) are given by:
\[ \psi_i = \sum_{j=1}^{2N} \tilde{\psi}_j \, (U_{\mathrm{proj}})_{j i} = \sum_{j=1}^{2N} \sum_{l=1}^{2N} \chi_l \, (\rho_{\mathrm{mat}})_{l j} \, (U_{\mathrm{proj}})_{j i} \]
where \( U_{\mathrm{proj}} = U V^{\mathrm{T}} \), and \( U \) and \( V \) are obtained from the singular value decomposition (SVD) of the overlap matrix \( S \):
\[ S = U \Sigma V^{\mathrm{T}} \]
Here, \( \rho_{\mathrm{mat}} \) contains the eigenvectors of the density matrix, and \( \{ \chi_i \} = \{ \phi_1, \dots, \phi_N, d_1, \dots, d_N \} \) is the complete basis.
The occupation associated with state \( \psi_i \) is given by \( |\psi|^2 \), where:
\[ \psi = \sum_{j=1}^{2N} \sum_{l=1}^{2N} \chi_l \, (\rho_{\mathrm{mat}})_{l j} \, \sqrt{\mathrm{occ}_j} \, (U_{\mathrm{proj}})_{j i} \]
Definition at line 756 of file dm_propagation.F90.
|
private |
Apply dissipation to the density matrix This subroutine applies the dissipation operator to the density matrix using Taylor expansion of the matrix exponential upto order 4.
\[ \overrightarrow{D}^{\mathrm{inter}}_{ab;cd} = \sum_{k=1}^{k_{\mathrm{max}}}\sum_{\substack{m,n = 0 \\ m \neq n}}^{N-1}\ & \gamma^{m \to n}_{k}\Bigg( \delta_{c,(m;k)}\ \delta_{(n;k),a} \ \delta_{b,(n;k)}\ \delta_{(m;k),d}- \frac{1}{2}\big[ & \delta_{c,a}\ \delta_{b,(m;k)}\ \Delta_{l_1}\ \delta_{(m;k),d}+ \delta_{b,d}\ \\ & \delta_{c,(m;k)}\ \Delta_{l_1}\delta_{(m;k),a}\big]\Bigg) \]
where \( \overrightarrow{D}^{\mathrm{inter}}\) is the interband dissipator using vectorized Lindblad formalism, \( \gamma^{m \to n}_{k}\) is the damping factor for the transition from band \(m\) to band \(n\) at $\f k[in] | dmp | Density matrix propagation instance |
[in] | bst | Adiabatic basis |
[in] | st | Time-dependent wfcs |
[in] | ik | K-point index |
[in] | dt | Time step |
[in,out] | nn | Density matrix |
Definition at line 867 of file dm_propagation.F90.