Kernels
This module provides kernels for the analytic continuation. In general, we have \(G_i = \sum_j K_{ij} H_j\). Note that for non-uniform \(\omega\)-grids, \(H = A \Delta\omega\).
At the moment, only the kernel for the continuation of \(G(\tau)\)
is implemented, i.e., the TauKernel
.
For any kernel the preblur version can be used by the PreblurKernel
(see Continuation of metallic solutions using preblur).
- class triqs_maxent.kernels.KernelSVD(K=None)[source]
Bases:
object
A kernel object with a singular value decomposition
- Parameters:
- Karray
the kernel matrix
- Attributes:
Methods
reduce_singular_space
([threshold])Reduce the singular space
svd
()Perform the SVD if not yet performed
- svd()[source]
Perform the SVD if not yet performed
Usually, this function does not need to be called by the user.
We have \(K = USV^\dagger\).
- property U
Get the matrix of left-singular vectors of the kernel
If not performed already, the SVD is done.
- property S
Get the vector of singular values of the kernel
If not performed already, the SVD is done.
- property V
Get the matrix of right-singular vectors of the kernel
If not performed already, the SVD is done.
- property K
The actual kernel matrix
- class triqs_maxent.kernels.Kernel[source]
Bases:
KernelSVD
The kernel for the analytic continuation
- Attributes:
K
The actual kernel matrix
K_delta
The kernel including a \(\Delta\omega\).
S
Get the vector of singular values of the kernel
U
Get the matrix of left-singular vectors of the kernel
V
Get the matrix of right-singular vectors of the kernel
- data_variable
Methods
Notify the kernel of a parameter change
reduce_singular_space
([threshold])Reduce the singular space
svd
()Perform the SVD if not yet performed
transform
(T_)multiply the kernel from the left with the matrix
T_
- property K_delta
The kernel including a \(\Delta\omega\).
Use this to get the reconstructed Green function as \(G_{rec} = K_{delta} A\). Note that it does not get rotated with
transform()
, i.e. it gives back the original \(G_{rec}\).
- class triqs_maxent.kernels.DataKernel(data_variable, omega, K)[source]
Bases:
Kernel
A kernel given by a matrix
- Parameters:
- data_variablearray
the array of the data variable, i.e. the variable that the input-data depends on; typically, one continues \(G(\tau)\), then this would correspond to the \(\tau\)-grid.
- omegaOmegaMesh
the frequency mesh
- Karray
the kernel matrix
- Attributes:
K
The actual kernel matrix
K_delta
The kernel including a \(\Delta\omega\).
S
Get the vector of singular values of the kernel
U
Get the matrix of left-singular vectors of the kernel
V
Get the matrix of right-singular vectors of the kernel
- data_variable
Methods
parameter_change
()Notify the kernel of a parameter change
reduce_singular_space
([threshold])Reduce the singular space
svd
()Perform the SVD if not yet performed
transform
(T_)multiply the kernel from the left with the matrix
T_
- class triqs_maxent.kernels.TauKernel(tau, omega, beta=None)[source]
Bases:
Kernel
A kernel for continuing \(G(\tau)\)
This kernel is defined as
\[K(\tau, \omega) = - \frac{\exp(-\tau\omega)}{1 + \exp(\beta \omega)}.\]With this, we have
\[G(\tau) = \int d\omega\, K(\tau, \omega) A(\omega).\]- Parameters:
- tauarray
the \(\tau\)-mesh where the data is given
- omegaOmegaMesh
the \(\omega\)-mesh where the spectral function should be calculated
- betafloat
the inverse temperature; if not given, it is taken to be the last
tau
-value
- Attributes:
K
The actual kernel matrix
K_delta
The kernel including a \(\Delta\omega\).
S
Get the vector of singular values of the kernel
U
Get the matrix of left-singular vectors of the kernel
V
Get the matrix of right-singular vectors of the kernel
data_variable
\(\tau\)
Methods
parameter_change
()Notify the kernel of a parameter change
reduce_singular_space
([threshold])Reduce the singular space
svd
()Perform the SVD if not yet performed
transform
(T_)multiply the kernel from the left with the matrix
T_
- property data_variable
\(\tau\)
- class triqs_maxent.kernels.IOmegaKernel(iomega, omega, beta=None)[source]
Bases:
Kernel
A kernel for continuing \(G(i\omega)\)
This kernel is defined as
\[K(i\omega, \omega) = \frac{1}{i\omega - \omega}\]With this, we have
\[G(i\omega) = \int d\omega\, K(i\omega, \omega) A(\omega).\]- Parameters:
- iomegaarray
the \(iomega\)-mesh where the data is given (as a real array, i.e. it is internally multiplied by
1.0j
)- omegaOmegaMesh
the \(\omega\)-mesh where the spectral function should be calculated
- betafloat
the inverse temperature; if not given, it is taken from the difference of the first two \(i\omega\) values
- Attributes:
K
The actual kernel matrix
K_delta
The kernel including a \(\Delta\omega\).
S
Get the vector of singular values of the kernel
U
Get the matrix of left-singular vectors of the kernel
V
Get the matrix of right-singular vectors of the kernel
data_variable
\(i\omega\)
Methods
parameter_change
()Notify the kernel of a parameter change
reduce_singular_space
([threshold])Reduce the singular space
svd
()Perform the SVD if not yet performed
transform
(T_)multiply the kernel from the left with the matrix
T_
- property data_variable
\(i\omega\)
- class triqs_maxent.kernels.PreblurKernel(K, b)[source]
Bases:
Kernel
A kernel for the preblur formalism
In the preblur formalism, the equation \(G = KA\) is replaced by \(G = KBH\), with a hidden image \(H\) and a blur matrix \(B\).
The dicretization of \(\omega\) has to be accounted for by including a \(\Delta\omega\). In fact, the total equation is \(G_i = \sum_{jk} K_{ij} \Delta\omega_j B_{jk} H_k\), where \(H_k\) already includes a \(\Delta\omega_k\).
Note that the
PreblurKernel
should always be used together with thePreblurA_of_H
.- Parameters:
- Attributes:
K
The actual kernel matrix
K_delta
The kernel including a \(\Delta\omega\).
S
Get the vector of singular values of the kernel
U
Get the matrix of left-singular vectors of the kernel
V
Get the matrix of right-singular vectors of the kernel
- data_variable
- omega
Methods
Notify the kernel of a parameter change
reduce_singular_space
([threshold])Reduce the singular space
svd
()Perform the SVD if not yet performed
transform
(T)multiply the kernel from the left with the matrix
T_
get_omega
set_omega