som
: Main Python module of SOM
som.Som
- class som.Som(rhs: Gf, errors: Gf, kind: str = 'FermionGf', norms=None, *, filtering_levels=None)[source]
Implementation of the Stochastic Optimization Method.
- __init__(rhs: Gf, errors: Gf, kind: str = 'FermionGf', norms=None, *, filtering_levels=None)[source]
- Parameters:
rhs (
triqs.gf.gf.Gf
) – Right hand side of the integral equation to be solved, defined ontriqs.gf.meshes.MeshImTime
,triqs.gf.meshes.MeshImFreq
ortriqs.gf.meshes.MeshLegendre
. The target shape ofrhs
must be \(M{\times}M\). If \(M>1\), only its diagonal matrix elements will be considered and used as input data for \(M\) independent continuation problems.errors – Either error bars \(\sigma_n\) (GF container of the same type and target shape as
rhs
) or a packed list of covariance matrices (Gf(mesh=MeshProduct(rhs.mesh, rhs.mesh), target_shape=[M])
with each target element corresponding to one covariance matrix).kind (
str
, optional) – Selection of the observable kind (and its respective integral kernel) to be used. One ofFermionGf
,FermionGfSymm
,BosonCorr
,BosonAutoCorr
,ZeroTemp
. Defaults toFermionGf
.norms (
float
orlist
[float
]) –Requested solution norms either as a list of \(M\) real numbers or as a single real number for all \(M\) continuation problems. For observable kinds
FermionGf
,FermionGfSymm
andZeroTemp
this parameter is optional and defaults to 1.0. ForBosonCorr
andBosonAutoCorr
it must be provided by the user.See also
When unknown a priori in the latter case, the norms can be approximately estimated from
rhs
usingestimate_boson_corr_spectrum_norms()
.filtering_levels (
float
orlist
[float
], optional) – Filtering levels for covariance matrices either as a list of \(M\) real numbers or as a single real number for all \(M\) continuation problems. Can only be specified when the covariance matrices are used. Defaults to 0.
- accumulate(**kwargs)
Accumulate particular solutions. This method can be called multiple times to accumulate the particular solutions incrementally.
- Parameters:
**kwargs – See below.
Main parameters
Parameter Name
Type
Default
Documentation
energy_window
(float,float)
–
Estimated lower and upper bounds of the spectrum.
max_time
int
-1 = infinite
Maximum runtime in seconds, use -1 to set infinite.
verbosity
int
2 on MPI rank 0, 0 otherwise.
Verbosity level (max level - 3).
t
int
50
Number of elementary updates per global update (\(T\)). Bigger values make the algorithm more ergodic.
cc_update
bool
False
Enable Consistent Constraints updates.
f
int
100
Number of global updates (\(F\)). Bigger values make the algorithm more ergodic.
l
int
2000
Number of particular solutions to accumulate (\(L\)); ignored if
adjust_l=True
. Bigger values reduce noise in the final solution / make it smoother.adjust_l
bool
False
Automatically adjust the number of particular solutions to accumulate (\(L\)).
make_histograms
bool
False
Accumulate histograms of \(\chi\), where \(\chi^2\) is the objective function.
Fine tuning options
Parameter Name
Type
Default
Documentation
random_seed
int
34788 + 928374 * MPI.rank
Seed for random number generator.
random_name
str
“”
Name of random number generator (MT19937 by default).
max_rects
int
60
Maximum number of rectangles in a particular solution (\(K_{max}\)).
min_rect_width
float
1e-3
Minimal width of a rectangle, in units of the energy window width.
min_rect_weight
float
1e-3
Minimal weight of a rectangle, in units of the requested solution norm.
t1
int
-1
Number of elementary updates in the first stage of a global update. When set to -1, the number of elementary updates will be chosen randomly for each global update from the \([1; T[\) range.
distrib_d_max
float
2.0
Maximal parameter of the power-law distribution function for the Metropolis algorithm.
gamma
float
2.0
Proposal probability parameter \(\gamma\).
cc_update_cycle_length
int
10
CC update: Number of proposed elementary updates between two successive CC updates (only during the first stage of a global update).
cc_update_max_iter
int
30
CC update: Maximum allowed number of height adjustment iterations.
cc_update_rect_norm_variation_tol
float
1e-3
CC update: The height adjustment procedure stops when variation of every rectangle’s norm between two consecutive iterations is below this value. This parameter is measured in units of the requested solution norm.
cc_update_height_penalty_max
float
1e3
CC update: Maximum value of the regularization parameters \(Q_0(k)\) that penalize negative heights. Measured in units of (energy window width) / (solution norm).
cc_update_height_penalty_divisor
float
10.0
CC update: Divisor used to reduce the regularization parameters \(Q_0(k)\) that penalize negative heights.
cc_update_der_penalty_init
float
1.0
CC update: Initial value of the regularization parameters \(Q_1(k)\) and \(Q_2(k)\) that penalize large derivatives of a solution. Measured in units of (energy window width)^2 / (solution norm) for \(Q_1(k)\) and in units of (energy window width)^3 / (solution norm) for \(Q_2(k)\).
cc_update_der_penalty_threshold
float
0.1
CC update: Sets the threshold value of the products \(|Q_1(k) A'(k)|\) and \(|Q_2(k) A''(k)|\), above which derivative regularization parameters \(Q_1(k)\) and \(Q_2(k)\) need to be reduced.
cc_update_der_penalty_increase_coeff
float
2.0
CC update: Coefficient used to increase the regularization parameters \(Q_1(k)\) and \(Q_2(k)\) that penalize large derivatives of a solution.
cc_update_der_penalty_limiter
float
1e3
CC update: Coefficient that limits growth of the regularization parameters \(Q_1(k)\) and \(Q_2(k)\) penalizing large derivatives of a solution.
adjust_l_range
(int,int)
(100,2000)
Automatic \(L\) adjustment: Search range for the number of particular solutions to accumulate.
adjust_l_good_chi
float
2.0
Automatic \(L\) adjustment: Maximal ratio \(\chi/\chi_\mathrm{min}\) for a particular solution to be considered good.
adjust_l_verygood_chi
float
4/3
Automatic \(L\) adjustment: Maximal ratio \(\chi/\chi_\mathrm{min}\) for a particular solution to be considered very good.
adjust_l_ratio
float
0.95
Automatic \(L\) adjustment: Critical ratio \(N_\mathrm{very good}/N_\mathrm{good}\) to stop the \(L\)-adjustment procedure.
hist_max
float
2.0
Right boundary of the histograms, in units of \(\chi_\mathrm{min}\) (left boundary is always set to \(\chi_\mathrm{min}\)).
hist_n_bins
int
100
Number of bins for the histograms.
- last_accumulate_parameters
Set of parameters used in the last call to
accumulate()
.
- accumulate_status
Status code of
accumulate()
on exit.0, if the accumulation has run until the end.
1, if it has run out of time.
2, if it has been stopped by receiving a signal.
- Type:
- adjust_f(**kwrags)
Automatically adjust the number of global updates \(F\).
This function scans a range of \(F\) (
f_range
), accumulatesl
particular solutions for each \(F\) and evaluates the following functional for each solution,\[\kappa = \frac{1}{N-1}\sum_{n=2}^N \frac{1}{2}\left[1 - \frac{\Re[\Delta(n)\Delta^*(n-1)]}{|\Delta(n)\Delta^*(n-1)|} \right],\]where \(\Delta(n)\) are discrepancies
\[\Delta(n) = \int_{\epsilon_\mathrm{min}}^{\epsilon_\mathrm{max}} K(n, \epsilon) A(\epsilon) d\epsilon - G(n).\]\(\kappa\) characterizes anti-correlation between adjacent discrepancies and can serve as a measure of ergodicity in Markov chains. If \(\kappa\) exceeds the value of parameter
kappa
for at least a half of accumulated solutions, \(F\) is considered sufficient.Deprecated since version 2.0: The automatic adjustment algorithm does not support continuation problems with full covariance matrices. If the
Som
object has been constructed using a covariance matrix,adjust_f()
will raise aRuntimeError
.- Parameters:
**kwargs – See below.
- Returns:
Selected value of \(F\).
- Return type:
Main parameters
Parameter Name
Type
Default
Documentation
energy_window
(float,float)
–
Estimated lower and upper bounds of the spectrum.
max_time
int
-1 = infinite
Maximum runtime in seconds, use -1 to set infinite.
verbosity
int
2 on MPI rank 0, 0 otherwise.
Verbosity level (max level - 3).
t
int
50
Number of elementary updates per global update (\(T\)). Bigger values make the algorithm more ergodic.
cc_update
bool
False
Enable Consistent Constraints updates.
f_range
(int,int)
(100,5000)
Search range for the number of global updates.
Fine tuning options
Parameter Name
Type
Default
Documentation
random_seed
int
34788 + 928374 * MPI.rank
Seed for random number generator.
random_name
str
“”
Name of random number generator (MT19937 by default).
max_rects
int
60
Maximum number of rectangles in a particular solution (\(K_{max}\)).
min_rect_width
float
1e-3
Minimal width of a rectangle, in units of the energy window width.
min_rect_weight
float
1e-3
Minimal weight of a rectangle, in units of the requested solution norm.
t1
int
-1
Number of elementary updates in the first stage of a global update. When set to -1, the number of elementary updates will be chosen randomly for each global update from the \([1; T[\) range.
distrib_d_max
float
2.0
Maximal parameter of the power-law distribution function for the Metropolis algorithm.
gamma
float
2.0
Proposal probability parameter \(\gamma\).
cc_update_cycle_length
int
10
CC update: Number of proposed elementary updates between two successive CC updates (only during the first stage of a global update).
cc_update_max_iter
int
30
CC update: Maximum allowed number of height adjustment iterations.
cc_update_rect_norm_variation_tol
float
1e-3
CC update: The height adjustment procedure stops when variation of every rectangle’s norm between two consecutive iterations is below this value. This parameter is measured in units of the requested solution norm.
cc_update_height_penalty_max
float
1e3
CC update: Maximum value of the regularization parameters \(Q_0(k)\) that penalize negative heights. Measured in units of (energy window width) / (solution norm).
cc_update_height_penalty_divisor
float
10.0
CC update: Divisor used to reduce the regularization parameters \(Q_0(k)\) that penalize negative heights.
cc_update_der_penalty_init
float
1.0
CC update: Initial value of the regularization parameters \(Q_1(k)\) and \(Q_2(k)\) that penalize large derivatives of a solution. Measured in units of (energy window width)^2 / (solution norm) for \(Q_1(k)\) and in units of (energy window width)^3 / (solution norm) for \(Q_2(k)\).
cc_update_der_penalty_threshold
float
0.1
CC update: Sets the threshold value of the products \(|Q_1(k) A'(k)|\) and \(|Q_2(k) A''(k)|\), above which derivative regularization parameters \(Q_1(k)\) and \(Q_2(k)\) need to be reduced.
cc_update_der_penalty_increase_coeff
float
2.0
CC update: Coefficient used to increase the regularization parameters \(Q_1(k)\) and \(Q_2(k)\) that penalize large derivatives of a solution.
cc_update_der_penalty_limiter
float
1e3
CC update: Coefficient that limits growth of the regularization parameters \(Q_1(k)\) and \(Q_2(k)\) penalizing large derivatives of a solution.
l
int
20
Number of particular solutions used to adjust \(F\).
kappa
float
0.25
Limiting value of \(\kappa\) used to adjust \(F\).
- compute_final_solution(good_chi_rel: float = 2.0, good_chi_abs: float = cmath.inf, verbosity: int = 0)
Signature : (float good_chi_rel = 2.0, float good_chi_abs = HUGE_VAL, int verbosity = 0) -> list[float]
Select particular solutions according to the standard SOM criterion and compute the final solution.
Selected solutions must satisfy \(\chi[A_j] \leq \chi_c^\mathrm{abs}\) and \(\chi[A_j] \leq \min_{j'}(\chi[A_{j'}]) \chi_c^\mathrm{rel}\).
- Parameters:
- Returns:
A length-\(M\) list of values of \(\chi^2\), one element per continuation problem.
- Return type:
- compute_final_solution_cc(**kwargs)
Select particular solutions and compute the final solution using the SOCC protocol.
Selected solutions must satisfy \(\chi[A_j] \leq \chi_c^\mathrm{abs}\) and \(\chi[A_j] \leq \min_{j'}(\chi[A_{j'}]) \chi_c^\mathrm{rel}\).
- Parameters:
**kwargs – See below.
- Returns:
A length-\(M\) list of values of \(\chi^2\), one element per continuation problem.
- Return type:
Main parameters
Parameter Name
Type
Default
Documentation
refreq_mesh
MeshReFreq or real 1D array_like
–
Grid of energy points \(\epsilon_k\) used in CC regularization procedure. Either a TRIQS real-frequency mesh object or a strictly ordered list of points.
verbosity
int
1 on MPI rank 0, 0 otherwise.
Verbosity level (max level - 2).
good_chi_rel
float
2.0
Maximal ratio \(\chi/\chi_\mathrm{min}\) for a particular solution to be selected. This criterion must be fulfilled together with the one set by
good_chi_abs
.good_chi_abs
float
infinity
Maximal value of \(\chi\) for a particular solution to be selected. This criterion must be fulfilled together with the one set by
good_chi_rel
.default_model
Real 1D array
[]
Optional default model of the spectral function evaluated at energy points of
refreq_mesh
(\(A_T(\epsilon_k)\)).default_model_weights
Real 1D array
[]
Weights determining how much deviations from
default_model
are penalized at each energy point (\(T_k\)).max_iter
int
50
Maximum allowed number of parameter adjustment iterations.
max_sum_abs_c
float
2.0
Stop parameter adjustment iterations when expansion coefficients \(c_j\) make the sum \(\sum_j |c_j|\) exceed this value.
Fine tuning options
Parameter Name
Type
Default
Documentation
ew_penalty_coeff
float
1.0
Coefficient of the term that penalizes large deviations from the equal-weight superposition.
amp_penalty_max
float
1e3
Maximum value of the regularization parameters that penalize negative values of the spectral function (\(\mathcal{Q}_k\)).
amp_penalty_divisor
float
10.0
Divisor used to reduce the regularization parameters that penalize negative values of the spectral function (\(\mathcal{Q}_k\)).
der_penalty_init
float
0.1
Initial value of the regularization parameters that penalize large derivatives of the solution (\(\mathcal{D}\)). Before this parameter is used, it is divided by the number of selected particular solutions.
der_penalty_coeff
float
2.0
Coefficient used to increase the regularization parameters that penalize large derivatives of the solution (\(f\)).
monitor
function
None
Monitor function called at each parameter adjustment iteration. It takes 4 arguments,
Current list of expansion coefficients \(c_j\);
Amplitudes of the spectrum and respective regularization parameters as a \((A_k, Q_k)\) pair;
Derivatives of the spectrum and respective regularization parameters as a \((A'_k, D_k)\) pair;
Second derivatives of the spectrum and respective regularization parameters as a \((A''_k, B_k)\) pair.
Returning
True
from the function stops the adjustment procedure.
- run(**kwargs)
Accumulate particular solutions and compute the final solution using the standard SOM criterion.
Calling this method is equivalent to calling
accumulate()
andcompute_final_solution()
withgood_chi_rel=adjust_l_good_chi
.Deprecated since version 2.0: Use
accumulate()
andcompute_final_solution()
instead.- Parameters:
**kwargs – See below.
Main parameters
Parameter Name
Type
Default
Documentation
energy_window
(float,float)
–
Estimated lower and upper bounds of the spectrum.
max_time
int
-1 = infinite
Maximum runtime in seconds, use -1 to set infinite.
verbosity
int
2 on MPI rank 0, 0 otherwise.
Verbosity level (max level - 3).
t
int
50
Number of elementary updates per global update (\(T\)). Bigger values make the algorithm more ergodic.
cc_update
bool
False
Enable Consistent Constraints updates.
f
int
100
Number of global updates (\(F\)). Bigger values make the algorithm more ergodic.
l
int
2000
Number of particular solutions to accumulate (\(L\)); ignored if
adjust_l=True
. Bigger values reduce noise in the final solution / make it smoother.adjust_l
bool
False
Automatically adjust the number of particular solutions to accumulate (\(L\)).
make_histograms
bool
False
Accumulate histograms of \(\chi\), where \(\chi^2\) is the objective function.
Fine tuning options
Parameter Name
Type
Default
Documentation
random_seed
int
34788 + 928374 * MPI.rank
Seed for random number generator.
random_name
str
“”
Name of random number generator (MT19937 by default).
max_rects
int
60
Maximum number of rectangles in a particular solution (\(K_{max}\)).
min_rect_width
float
1e-3
Minimal width of a rectangle, in units of the energy window width.
min_rect_weight
float
1e-3
Minimal weight of a rectangle, in units of the requested solution norm.
t1
int
-1
Number of elementary updates in the first stage of a global update. When set to -1, the number of elementary updates will be chosen randomly for each global update from the \([1; T[\) range.
distrib_d_max
float
2.0
Maximal parameter of the power-law distribution function for the Metropolis algorithm.
gamma
float
2.0
Proposal probability parameter \(\gamma\).
cc_update_cycle_length
int
10
CC update: Number of proposed elementary updates between two successive CC updates (only during the first stage of a global update).
cc_update_max_iter
int
30
CC update: Maximum allowed number of height adjustment iterations.
cc_update_rect_norm_variation_tol
float
1e-3
CC update: The height adjustment procedure stops when variation of every rectangle’s norm between two consecutive iterations is below this value. This parameter is measured in units of the requested solution norm.
cc_update_height_penalty_max
float
1e3
CC update: Maximum value of the regularization parameters \(Q_0(k)\) that penalize negative heights. Measured in units of (energy window width) / (solution norm).
cc_update_height_penalty_divisor
float
10.0
CC update: Divisor used to reduce the regularization parameters \(Q_0(k)\) that penalize negative heights.
cc_update_der_penalty_init
float
1.0
CC update: Initial value of the regularization parameters \(Q_1(k)\) and \(Q_2(k)\) that penalize large derivatives of a solution. Measured in units of (energy window width)^2 / (solution norm) for \(Q_1(k)\) and in units of (energy window width)^3 / (solution norm) for \(Q_2(k)\).
cc_update_der_penalty_threshold
float
0.1
CC update: Sets the threshold value of the products \(|Q_1(k) A'(k)|\) and \(|Q_2(k) A''(k)|\), above which derivative regularization parameters \(Q_1(k)\) and \(Q_2(k)\) need to be reduced.
cc_update_der_penalty_increase_coeff
float
2.0
CC update: Coefficient used to increase the regularization parameters \(Q_1(k)\) and \(Q_2(k)\) that penalize large derivatives of a solution.
cc_update_der_penalty_limiter
float
1e3
CC update: Coefficient that limits growth of the regularization parameters \(Q_1(k)\) and \(Q_2(k)\) penalizing large derivatives of a solution.
adjust_l_range
(int,int)
(100,2000)
Automatic \(L\) adjustment: Search range for the number of particular solutions to accumulate.
adjust_l_good_chi
float
2.0
Automatic \(L\) adjustment: Maximal ratio \(\chi/\chi_\mathrm{min}\) for a particular solution to be considered good.
adjust_l_verygood_chi
float
4/3
Automatic \(L\) adjustment: Maximal ratio \(\chi/\chi_\mathrm{min}\) for a particular solution to be considered very good.
adjust_l_ratio
float
0.95
Automatic \(L\) adjustment: Critical ratio \(N_\mathrm{very good}/N_\mathrm{good}\) to stop the \(L\)-adjustment procedure.
hist_max
float
2.0
Right boundary of the histograms, in units of \(\chi_\mathrm{min}\) (left boundary is always set to \(\chi_\mathrm{min}\)).
hist_n_bins
int
100
Number of bins for the histograms.
- clear()
Signature : () -> None
Discard all accumulated particular solutions, histograms and final solutions.
- observable_kind
Kind of the physical observable selected upon construction, one of
FermionGf
,FermionGfSymm
,BosonCorr
,BosonAutoCorr
,ZeroTemp
.- Type:
- particular_solutions(i: int)
Signature : (int i) -> list[std::pair<configuration, double>]
Particular solutions and their respective values of the objective function \(\chi^2\) for the i-th diagonal element of the observable accumulated on this MPI rank.
- Parameters:
i (
int
) – Index of the diagonal matrix element of the observable.- Return type:
list
[tuple
[Configuration
,float
)]]
- objf_min
Minimum of the objective function \(\chi^2\) over all accumulated particular solutions (one value per diagonal matrix element of the observable).
- Type:
- solution(i: int)
Signature : (int i) -> Configuration
Final solution for the i-th diagonal element of the observable.
- Parameters:
i (
int
) – Index of the diagonal matrix element of the observable.- Return type:
- solutions
Final solutions, one per diagonal element of the observable.
- Type:
- objf(i: int)
Signature : (int i) -> float
Value of the objective function \(\chi^2\) of the final solution for the i-th diagonal element of the observable.
- objf_list
Values of the objective function \(\chi^2\) of the final solutions, one value per diagonal element of the observable.
- histogram(i: int)
Signature : (int i) -> std::optional<histogram>
Accumulated \(\chi\) histogram for the i-th diagonal element of the observable. Returns
None
if histograms have not been accumulated.- Parameters:
i (
int
) – Index of the diagonal matrix element of the observable.- Return type:
Free functions
- som.fill_refreq()
Signature : (triqs::gfs::gf_view<refreq> g_w, SomCore cont, bool with_binning = true) -> None Fill a Green’s function object with the real-frequency version of the continued observable from a computed SOM solution.
Parameters:
- G_w:
triqs.gf.gf.Gf
defined ontriqs.gf.meshes.MeshReFreq
, target Green’s function object.- Cont:
Som
, Analytic continuation object.- With_binning:
bool
, Use binning while projecting onto the real-frequency mesh.
- som.compute_tail()
Signature : (int max_order, SomCore cont) -> nda::array<dcomplex, 3> Extract high-frequency expansion coefficients from a computed SOM solution.
Parameters:
Returns: 3D complex NumPy array of shape
(max_order+1, cont.dim, cont.dim)
with the high-frequency expansion coefficients.
- som.reconstruct()
Signature : (triqs::gfs::gf_view<imtime> g, SomCore cont) -> None Reconstruct a physical observable in the imaginary-time representation from a computed SOM solution.
Parameters:
- G:
triqs.gf.gf.Gf
defined ontriqs.gf.meshes.MeshImTime
, target Green’s function object.- Cont:
Som
, Analytic continuation object.
- som.estimate_boson_corr_spectrum_norms(chi: Gf) List[float] [source]
Given a correlator \(\chi\) of bosonic or boson-like (fermion-number-conserving) operators, estimates the corresponding spectrum normalization constants \(\mathcal{N}\).
Depending on the mesh \(\chi\) is defined on, one of the following expressions is used.
Imaginary frequencies: \(\mathcal{N} = \pi\chi(i\Omega=0)\).
Imaginary time: \(\mathcal{N} = \pi\int_0^\beta d\tau \chi(\tau)\).
Legendre polynomial basis coefficients: \(\mathcal{N} = \pi\chi(\ell=0)\).
- Parameters:
chi (
triqs.gf.gf.Gf
) – The correlator \(\chi\).- Returns:
A list of estimated spectrum normalization constants, one constant per diagonal matrix element of \(\chi\).
- Return type:
- som.count_good_solutions(hist: Histogram, good_chi_rel: float = 2.0, good_chi_abs: float = inf) int [source]
Given a histogram of values \(\chi\) for the objective function \(\chi^2\), counts the number of such solutions that \(\chi \leq\)
good_chi_abs
and \(\chi/\chi_\mathrm{min} \leq\)good_chi_rel
.- Parameters:
hist (
triqs.stat.histograms.Histogram
) – Histogram to analyze.good_chi_rel (
float
, optional) – \(\chi/\chi_\mathrm{min}\) threshold for good solutions.good_chi_abs (
float
, optional) – \(\chi\) threshold for good solutions.
- Return type:
Auxiliary types
- class som.Rectangle
The rectangular function of energy, basis element of a spectral function.
- center
float: Center of the rectangle \(c\).
- width
float: Width of the rectangle \(w\).
- height
float: Height of the rectangle \(h\).
- norm
float: Norm (area) of the rectangle \(hw\).