dyson.solvers.dynamic.corrvec#
Correction vector Green’s function solver [1].
Classes
|
Correction vector Green's function solver. |
- class dyson.solvers.dynamic.corrvec.CorrectionVector(*args: Any, **kwargs: Any)[source]#
Bases:
DynamicSolverCorrection vector Green’s function solver.
- Parameters:
matvec – The matrix-vector operation for the self-energy supermatrix.
diagonal – The diagonal of the self-energy supermatrix.
nphys – The number of physical degrees of freedom.
grid – Real frequency grid upon which to evaluate the Green’s function.
- classmethod from_self_energy(static: Array, self_energy: Lehmann, overlap: Array | None = None, **kwargs: Any) CorrectionVector[source]#
Create a solver from a self-energy.
- Parameters:
static – Static part of the self-energy.
self_energy – Self-energy.
overlap – Overlap matrix for the physical space.
kwargs – Additional keyword arguments for the solver.
- Returns:
Solver instance.
- classmethod from_expression(expression: BaseExpression, **kwargs: Any) CorrectionVector[source]#
Create a solver from an expression.
- Parameters:
expression – Expression to be solved.
kwargs – Additional keyword arguments for the solver.
- Returns:
Solver instance.
- matvec_dynamic(vector: Array, grid_index: int) Array[source]#
Perform the matrix-vector operation for the dynamic self-energy supermatrix.
\[\mathbf{x}_\omega = \left(\omega - \mathbf{H} - i\eta\right) \mathbf{r}\]- Parameters:
vector – The vector to operate on.
grid_index – Index of the real frequency grid.
- Returns:
The result of the matrix-vector operation.
- matdiv_dynamic(vector: Array, grid_index: int) Array[source]#
Approximately perform a matrix-vector division for the dynamic self-energy supermatrix.
\[\mathbf{x}_\omega = \frac{\mathbf{r}}{\omega - \mathbf{H} - i\eta}\]- Parameters:
vector – The vector to operate on.
grid_index – Index of the real frequency grid.
- Returns:
The result of the matrix-vector division.
Notes
The inversion is approximated using the diagonal of the matrix.
- get_state_bra(orbital: int) Array[source]#
Get the bra vector corresponding to a fermion operator acting on the ground state.
- Parameters:
orbital – Orbital index.
- Returns:
Bra vector.
- get_state_ket(orbital: int) Array[source]#
Get the ket vector corresponding to a fermion operator acting on the ground state.
- Parameters:
orbital – Orbital index.
- Returns:
Ket vector.
- kernel(*args: Any, **kwargs: Any) Any#
Run the solver.
- property matvec: Callable[[Array], Array]#
Get the matrix-vector operation for the self-energy supermatrix.
- property diagonal: Array#
Get the diagonal of the self-energy supermatrix.
- property grid: RealFrequencyGrid#
Get the real frequency grid.