TRIQS/triqs_modest 3.3.0
Brillouin zone summation
|
Namespaces | |
namespace | dft_code |
Classes | |
struct | atomic_orbs |
Info on an atomic shell. More... | |
class | band_dispersion |
The one-body dispersion as a function of momentum. More... | |
class | checkpoint |
A checkpoint manager for logging data after each DMFT iteration. More... | |
class | dc_solver |
Double counting "solver" implements the double counting correction for DFT+DMFT, which is a phenomenlogical introduced double counting the interactions already taken into account at the mean-field level within DFT. This class implements several double counting formulas (all of which are functions of the density) relevant for different scenarios. More... | |
class | downfolding_projector |
The projector that downfolds the one-body dispersion (ν) onto local orbitals (m). More... | |
class | embedding |
The embedding class. More... | |
class | ibz_symmetry_ops |
ibz symmetry operations More... | |
struct | initial_data |
struct | iteration_data |
class | local_space |
Describe the atomic orbitals within downfolded \(\mathcal{C}\) space. More... | |
struct | one_body_elements_on_grid |
A one-body elements struct where all of the underlying data exists on a fixed momentum grid. More... | |
struct | one_body_elements_tb |
struct | spectral_function_kw |
Returns Tr (A) [σ,k,ω] for all k points in obe grid and all omega in Sigma mesh. More... | |
struct | spectral_function_w |
Typedefs | |
using | block_gf_imfreq_t = block_gf< imfreq, matrix_valued > |
using | block_gf_imtime_t = block_gf< imtime, matrix_valued > |
using | block_mat_t = std::vector< nda::matrix< dcomplex > > |
Enumerations | |
enum class | DFTCode { Wien2k , QuantumEspresso , VASP , Elk , W90 , Hk } |
enum class | ReadMode { Correlated , ThetaProjectors , Bands } |
Enumerate the different reading modes for the obe factory functions. More... | |
enum class | spin_kind_e { Polarized , NonPolarized , NonColinear } |
Kind of σ index. More... | |
Functions | |
std::vector< std::vector< long > > | analyze_degenerate_blocks (block_gf< imfreq, matrix_valued > const &Gimp, double threshold=1.e-5) |
Find the generate blocks of a BlockGf by analyzing G(τ=0) or G(iω₀) using the union-find algorithm. | |
std::pair< double, double > | dc_formulas (std::string const method, double const N_tot, double const N_sigma, long const n_orb, double const U, double const J) |
double counting formulas parameterized by density, U, and J | |
template<typename Mesh > | |
double | density (one_body_elements_on_grid const &obe, double mu, block2_gf< Mesh, matrix_valued > const &Sigma_dynamic, nda::array< nda::matrix< dcomplex >, 2 > const &Sigma_static) |
Compute the density of the lattice Green's function with a self-energy using Woodbury. | |
template<typename Mesh > | |
double | density (one_body_elements_tb const &obe, double mu, block2_gf< Mesh, matrix_valued > const &Sigma_dynamic, nda::array< nda::matrix< dcomplex >, 2 > const &Sigma_static, triqs::lattice::bz_int_options const &opt) |
Compute the density of the lattice Green's function with a self-energy. | |
double | density_nk (one_body_elements_on_grid const &obe, double mu, double beta) |
Compute number of particles n = ∑f(β(e(k)-μ)) | |
double | density_nk_matrix_valued_impl (one_body_elements_on_grid const &obe, double mu, double beta) |
std::pair< nda::array< std::vector< long >, 2 >, nda::array< nda::matrix< dcomplex >, 2 > > | discover_symmetries (nda::array< nda::matrix< dcomplex >, 2 > const &Hloc0, std::vector< atomic_orbs > const &atomic_shells, double block_threshold, bool diagonalize_hloc) |
Find symmetries of the R = 0 component of a Hamiltonian to determine a GF block structure. | |
std::pair< nda::array< nda::matrix< double >, 2 >, nda::matrix< double > > | double_counting (nda::array< nda::matrix< dcomplex >, 2 > const &density_matrix, double U_int, double J_hund, std::string const method) |
compute double counting correction for a dc_type (method) from the density matrix of a Green's function. | |
nda::matrix< double > | double_counting_energy_dc (nda::array< nda::matrix< dcomplex >, 2 > const &density_matrix, double U_int, double J_hund, std::string const method) |
nda::array< nda::matrix< double >, 2 > | double_counting_sigma_dc (nda::array< nda::matrix< dcomplex >, 2 > const &density_matrix, double U_int, double J_hund, std::string const method) |
template<typename Mesh > | |
double | find_chemical_potential (double const target_density, one_body_elements_on_grid const &obe, block2_gf< Mesh, matrix_valued > const &Sigma_dynamic, nda::array< nda::matrix< dcomplex >, 2 > const &Sigma_static, std::string method="dichotomy", double precision=1.e-5, bool verbosity=true) |
Find the chemical potenital from the local Green's function and self-energy given a target density. | |
double | find_chemical_potential (double const target_density, one_body_elements_on_grid const &obe, double beta, std::string method="dichotomy", double precision=1.e-5, bool verbosity=true) |
Find the chemical potenital from the local Green's function given a target density. | |
template<typename Mesh > | |
double | find_chemical_potential (double const target_density, one_body_elements_tb const &obe, block2_gf< Mesh, matrix_valued > const &Sigma_dynamic, nda::array< nda::matrix< dcomplex >, 2 > const &Sigma_static, triqs::lattice::bz_int_options const &opt, std::string method="dichotomy", double precision=1.e-5, bool verbosity=true) |
Find the chemical potenital from the local Green's function and self-energy given a target density. | |
template<typename Mesh > | |
double | find_chemical_potential (double const target_density, one_body_elements_tb const &obe, Mesh const &mesh, triqs::lattice::bz_int_options const &opt, std::string method="dichotomy", double precision=1.e-5, bool verbosity=true) |
Find the chemical potenital from the local Green's function and self-energy given a target density. | |
one_body_elements_tb | fold (lattice::superlattice const &sl, one_body_elements_tb const &obe) |
one_body_elements_tb | fold (superlattice const &sl, one_body_elements_tb const &obe) |
Convert a tight binding Hamiltonian to its superlattice equivalent. | |
auto | format_as (embedding::imp_block_t const &p) |
std::pair< double, nda::vector< double > > | get_total_density (nda::array< nda::matrix< dcomplex >, 2 > const &density_matrix) |
template block2_gf< dlr_imfreq, matrix_valued > | gloc (dlr_imfreq const &mesh, one_body_elements_on_grid const &obe, double mu) |
template block2_gf< imfreq, matrix_valued > | gloc (imfreq const &mesh, one_body_elements_on_grid const &obe, double mu) |
template<typename Mesh > | |
block2_gf< Mesh, matrix_valued > | gloc (Mesh const &mesh, one_body_elements_tb const &obe, double mu, triqs::lattice::bz_int_options const &opt) |
Compute the local Green's function without a self-energy. | |
template block2_gf< dlr_imfreq, matrix_valued > | gloc (one_body_elements_on_grid const &one_body, double mu, block2_gf< dlr_imfreq, matrix_valued > const &Sigma_dynamic, nda::array< nda::matrix< dcomplex >, 2 > const &Sigma_static) |
template block2_gf< imfreq, matrix_valued > | gloc (one_body_elements_on_grid const &one_body, double mu, block2_gf< imfreq, matrix_valued > const &Sigma_dynamic, nda::array< nda::matrix< dcomplex >, 2 > const &Sigma_static) |
template<typename Mesh > | |
block2_gf< Mesh, matrix_valued > | gloc (one_body_elements_tb const &obe, double mu, block2_gf< Mesh, matrix_valued > const &Sigma_dynamic, nda::array< nda::matrix< dcomplex >, 2 > const &Sigma_static, triqs::lattice::bz_int_options const &opt) |
Compute the local Green's function without a self-energy. | |
void | h5_read (h5::group g, std::string const &name, atomic_orbs &shell) |
void | h5_read (h5::group g, std::string const &name, band_dispersion &bd) |
void | h5_read (h5::group g, std::string const &name, downfolding_projector &proj) |
void | h5_read (h5::group g, std::string const &name, embedding &embed) |
void | h5_read (h5::group g, std::string const &name, initial_data &data) |
void | h5_read (h5::group g, std::string const &name, iteration_data &data) |
void | h5_read (h5::group g, std::string const &name, local_space &ls) |
void | h5_read (h5::group g, std::string const &name, one_body_elements_on_grid &x) |
void | h5_read (h5::group g, std::string const &name, spin_kind_e &spin_kind) |
void | h5_write (h5::group g, std::string const &name, atomic_orbs const &shell) |
void | h5_write (h5::group g, std::string const &name, band_dispersion const &bd) |
void | h5_write (h5::group g, std::string const &name, downfolding_projector const &proj) |
void | h5_write (h5::group g, std::string const &name, embedding const &embed) |
void | h5_write (h5::group g, std::string const &name, initial_data const &data) |
void | h5_write (h5::group g, std::string const &name, iteration_data const &data) |
void | h5_write (h5::group g, std::string const &name, local_space const &ls) |
void | h5_write (h5::group g, std::string const &name, one_body_elements_on_grid const &x) |
void | h5_write (h5::group g, std::string const &name, spin_kind_e const &spin_kind) |
nda::array< nda::matrix< dcomplex >, 2 > | Hloc (std::vector< tb_hamiltonian > const &H_sigma, std::vector< atomic_orbs > const &atomic_shells) |
Compute Hloc = H(R=0) given n_sigma tight_binding Hamiltonians. | |
template block_gf< dlr_imfreq, matrix_valued > | hybridization (std::vector< nda::matrix< dcomplex > > const &epsilon_levels, block_gf< dlr_imfreq, matrix_valued > const &Gloc) |
template block_gf< dlr_imfreq, matrix_valued > | hybridization (std::vector< nda::matrix< dcomplex > > const &epsilon_levels, block_gf< dlr_imfreq, matrix_valued > const &Gloc, block_gf< dlr_imfreq, matrix_valued > const &Sigma_dynamic, std::vector< nda::matrix< dcomplex > > const &Sigma_static) |
template block_gf< imfreq, matrix_valued > | hybridization (std::vector< nda::matrix< dcomplex > > const &epsilon_levels, block_gf< imfreq, matrix_valued > const &Gloc) |
template block_gf< imfreq, matrix_valued > | hybridization (std::vector< nda::matrix< dcomplex > > const &epsilon_levels, block_gf< imfreq, matrix_valued > const &Gloc, block_gf< imfreq, matrix_valued > const &Sigma_dynamic, std::vector< nda::matrix< dcomplex > > const &Sigma_static) |
template<typename Mesh > | |
block_gf< Mesh, matrix_valued > | hybridization (std::vector< nda::matrix< dcomplex > > const &epsilon_levels, block_gf< Mesh, matrix_valued > const &Gloc) |
Compute the hybridization function from the effective impurity levels and the local Green's function. | |
template<typename Mesh > | |
block_gf< Mesh, matrix_valued > | hybridization (std::vector< nda::matrix< dcomplex > > const &epsilon_levels, block_gf< Mesh, matrix_valued > const &Gloc, block_gf< Mesh, matrix_valued > const &Sigma_dynamic, std::vector< nda::matrix< dcomplex > > const &Sigma_static) |
Compute the hybridization function from the effective impurity levels, the local Green's function, and the impurity self-energy. | |
nda::array< nda::matrix< dcomplex >, 2 > | impurity_levels (one_body_elements_on_grid const &obe) |
Compute the local impurity levels from the single-particle dispersion. | |
nda::array< nda::matrix< dcomplex >, 2 > | impurity_levels (one_body_elements_tb const &obe) |
Compute the atomic (impurity) levels from an obe. | |
INSTANTIATE (dlr_imfreq) | |
INSTANTIATE (imfreq) | |
INSTANTIATE (refreq) | |
constexpr auto | lattice_gf (auto const &mesh, auto const &obe, auto const &Proj, auto const &mu, auto const &embedding_decomp, auto const &Sigma_w, auto const &broadening) |
nda::array< dcomplex, 4 > | load_rotate_and_format_projectors (std::string const &filename, ReadMode mode, std::vector< cmat_t > const &rot_mats, std::vector< long > const &atom_decomp) |
Read projectors using ReadMode, rotate to local frame, and embed them in the M space. (internal) | |
embedding | make_embedding_impl (local_space const &C_space, nda::array< std::vector< long >, 2 > const &block_decomposition, std::optional< std::vector< long > > const &atom_to_imp) |
std::pair< one_body_elements_on_grid, embedding > | make_embedding_with_clusters (one_body_elements_on_grid obe, std::vector< std::vector< long > > const &atom_partition) |
embedding | make_embedding_with_equivalences (local_space const &C_space) |
embedding | make_embedding_with_no_equivalences (local_space const &C_space) |
one_body_elements_tb | make_obe_from_tb (std::vector< tb_hamiltonian > const tb_H_sigma, spin_kind_e spin_kind, std::vector< atomic_orbs > atomic_shells) |
Helper to contruct and return an OBE_tb object given a list of tb_Hamiltonians of length n_sigma. | |
one_body_elements_tb | one_body_elements_from_wannier90 (std::string const &wannier_file_path, spin_kind_e spin_kind, std::vector< atomic_orbs > atomic_shells) |
Construct a obe_tb from Wannier90 in the case of a single spin index. | |
one_body_elements_tb | one_body_elements_from_wannier90 (std::string const &wannier_file_path_up, std::string const &wannier_file_path_dn, spin_kind_e spin_kind, std::vector< atomic_orbs > atomic_shells) |
Construct a obe_tb from Wannier90 in the case with separate spin up/spin down channels. | |
std::ostream & | operator<< (std::ostream &out, band_dispersion const &x) |
std::ostream & | operator<< (std::ostream &out, downfolding_projector const &x) |
std::ostream & | operator<< (std::ostream &out, embedding const &E) |
std::ostream & | operator<< (std::ostream &out, ibz_symmetry_ops const &ibz) |
std::ostream & | operator<< (std::ostream &out, local_space const &bd) |
std::ostream & | operator<< (std::ostream &out, one_body_elements_on_grid const &) |
one_body_elements_on_grid | permute_local_space (std::vector< std::vector< long > > const &atom_partition, one_body_elements_on_grid const &obe) |
spectral_function_w | projected_spectral_function (one_body_elements_on_grid const &obe_theta, downfolding_projector const &Proj, double mu, block2_gf< mesh::refreq, matrix_valued > const &Sigma_w, double broadening=0.01) |
Compute the atom- and orbital-resolved spectral function (interacting density of states). | |
std::vector< atomic_orbs > | read_atomic_shells (auto const &filename, ReadMode mode) |
Read atomic shells according to ReadMode. (internal) | |
std::tuple< nda::array< dcomplex, 4 >, nda::matrix< long >, nda::array< double, 1 > > | read_bands_and_weights (std::string filename, ReadMode mode) |
Read band dispersion and k-weights according to ReadMode. (internal) | |
ibz_symmetry_ops | read_ibz_symmetry_ops (auto const &filename, ReadMode mode) |
Construct the ibz_symmetry_ops according to ReadMode. | |
std::vector< cmat_t > | read_rotation_matrices (std::string const &filename, ReadMode mode) |
Read rotation matrices from hdf5. (internal) | |
nda::array< nda::matrix< dcomplex >, 1 > | read_spherical_to_dft_basis (std::string dft, std::vector< atomic_orbs > const &atomic_shells) |
Prepare the spherical Ylm to DFT orbital basis rotations. (internal) | |
spin_kind_e | read_spin_kind (auto const &filename) |
Setup spin_kind enum. (internal) | |
ibz_symmetry_ops | rotate_local_basis (nda::array< nda::matrix< dcomplex >, 2 > const &U, ibz_symmetry_ops const &x) |
Change basis. | |
one_body_elements_on_grid | rotate_local_basis (nda::array< nda::matrix< dcomplex >, 2 > const &U, one_body_elements_on_grid const &x) |
Rotates the local basis of the downfolding projector. | |
long | sigma_to_data_idx (spin_kind_e spin_kind, long sigma) |
spectral_function_kw | spectral_function_on_high_symmetry_path (one_body_elements_on_grid const &obe, double mu, block2_gf< mesh::refreq, matrix_valued > const &Sigma_w, double broadening=0.01) |
Compute momentum-resolved spectral function A(k, ω) along high-symmetry path. | |
nda::matrix< dcomplex > | svd (const nda::matrix< dcomplex > &A) |
utility to flatten a nested vector (move to utils/ ?) | |
block_gf< imfreq, matrix_valued > | symmetrize_gf (block_gf< imfreq, matrix_valued > const &Gin, std::vector< std::vector< long > > degenerate_blocks) |
Symmetrize the blocks of a Green's function given a list of it's degenerate blocks. | |
Embedding factories functions | |
Factory functions to create the embedding class for different embedding scenarios. Typically, one will create the embedding from the local space. | |
embedding | make_embedding (local_space const &C_space, bool use_atom_equivalences=true) |
Make an embedding from the local space. | |
Local Green's function using a fixed k-grid | |
template<typename Mesh > | |
block2_gf< Mesh, matrix_valued > | gloc (one_body_elements_on_grid const &obe, double mu, block2_gf< Mesh, matrix_valued > const &Sigma_dynamic, nda::array< nda::matrix< dcomplex >, 2 > const &Sigma_static) |
compute G𝓒 local Green's function on Mesh(MxM) | |
template<typename Mesh > | |
block2_gf< Mesh, matrix_valued > | gloc (Mesh const &mesh, one_body_elements_on_grid const &obe, double mu) |
Compute the local Green's function without a self-energy. | |
OBE factories using a fixed grid | |
Factory functions to create one_body_elements_on_grid | |
std::pair< double, one_body_elements_on_grid > | one_body_elements_from_dft_converter (std::string const &filename, double threshold, bool diagonalize_hloc) |
Prepare one-body elements for a DMFT calculation. | |
one_body_elements_on_grid | one_body_elements_with_theta_projectors (std::string const &filename, one_body_elements_on_grid const &obe) |
Prepare one-body elements with the Θ projectors. | |
one_body_elements_on_grid | one_body_elements_on_high_symmetry_path (std::string const &filename, one_body_elements_on_grid const &obe) |
Prepare one-body elements along high-symmetry k-path. | |
using triqs::modest::block_gf_imfreq_t = typedef block_gf<imfreq, matrix_valued> |
Definition at line 120 of file checkpoint.hpp.
using triqs::modest::block_gf_imtime_t = typedef block_gf<imtime, matrix_valued> |
Definition at line 121 of file checkpoint.hpp.
using triqs::modest::block_mat_t = typedef std::vector<nda::matrix<dcomplex> > |
Definition at line 122 of file checkpoint.hpp.
|
strong |
Enumerator | |
---|---|
Wien2k | |
QuantumEspresso | |
VASP | |
Elk | |
W90 | |
Hk |
Definition at line 10 of file dft_code_specific.hpp.
|
strong |
Enumerate the different reading modes for the obe factory functions.
Enumerator | |
---|---|
Correlated | |
ThetaProjectors | |
Bands |
Definition at line 83 of file loaders.cpp.
|
strong |
Kind of σ index.
Enumerator | |
---|---|
Polarized | |
NonPolarized | |
NonColinear |
Definition at line 17 of file local_space.hpp.
double triqs::modest::density_nk | ( | one_body_elements_on_grid const & | obe, |
double | mu, | ||
double | beta | ||
) |
Compute number of particles n = ∑f(β(e(k)-μ))
Definition at line 29 of file density.cpp.
double triqs::modest::density_nk_matrix_valued_impl | ( | one_body_elements_on_grid const & | obe, |
double | mu, | ||
double | beta | ||
) |
Definition at line 10 of file density.cpp.
std::pair< nda::array< std::vector< long >, 2 >, nda::array< nda::matrix< dcomplex >, 2 > > triqs::modest::discover_symmetries | ( | nda::array< nda::matrix< dcomplex >, 2 > const & | Hloc0, |
std::vector< atomic_orbs > const & | atomic_shells, | ||
double | block_threshold, | ||
bool | diagonalize_hloc | ||
) |
Find symmetries of the R = 0 component of a Hamiltonian to determine a GF block structure.
Disovers (approximate) irreducible symmetries for Green's function from the non-interacting part of the local Hamiltonian (H0 = ∑k P(k) Hνν' P†(k) ), which represents the block structure of the TRIQS Gf.
Discovers (approximate) irreducible symmetries for Green's function from the non-interacting part of the local Hamiltonian (H0 = ∑k P(k) Hνν' P†(k) ), which represents the block structure of the TRIQS Gf.
Hloc0 | the R = 0 part of the Hamiltonian as a vector of [n_atoms, n_sigma] |
atomic_shells | the list of atomic shells used to index Hloc |
block_threshold | the threshold of accuracy at which a symmetry is considered found |
diagonalize_hloc | whether or not to diagonalize hloc |
Definition at line 54 of file loaders.cpp.
std::pair< nda::array< nda::matrix< double >, 2 >, nda::matrix< double > > triqs::modest::double_counting | ( | nda::array< nda::matrix< dcomplex >, 2 > const & | density_matrix, |
double | U_int, | ||
double | J_hund, | ||
std::string const | method | ||
) |
compute double counting correction for a dc_type (method) from the density matrix of a Green's function.
density_matrix | the density matrix [orbital, spin] indices |
U_int | Coulomb interaction parameter |
J_hund | Hund's coupling interaction parameter |
method | dc_formula (sFLL, cFLL, sAMF, cAMF, cHeld) |
Definition at line 51 of file double_counting.cpp.
nda::matrix< double > triqs::modest::double_counting_energy_dc | ( | nda::array< nda::matrix< dcomplex >, 2 > const & | density_matrix, |
double | U_int, | ||
double | J_hund, | ||
std::string const | method | ||
) |
Definition at line 128 of file double_counting.cpp.
nda::array< nda::matrix< double >, 2 > triqs::modest::double_counting_sigma_dc | ( | nda::array< nda::matrix< dcomplex >, 2 > const & | density_matrix, |
double | U_int, | ||
double | J_hund, | ||
std::string const | method | ||
) |
Definition at line 108 of file double_counting.cpp.
one_body_elements_tb triqs::modest::fold | ( | lattice::superlattice const & | sl, |
one_body_elements_tb const & | obe | ||
) |
Definition at line 145 of file obe_tb.cpp.
one_body_elements_tb triqs::modest::fold | ( | superlattice const & | sl, |
one_body_elements_tb const & | obe | ||
) |
Convert a tight binding Hamiltonian to its superlattice equivalent.
sl | The superlattice object containing its lattice vectors and locations of cluster points. |
obe | A one_body_elements object containing the tb_hamiltonian. |
auto triqs::modest::format_as | ( | embedding::imp_block_t const & | p | ) |
Definition at line 91 of file embedding.cpp.
std::pair< double, nda::vector< double > > triqs::modest::get_total_density | ( | nda::array< nda::matrix< dcomplex >, 2 > const & | density_matrix | ) |
Definition at line 90 of file double_counting.cpp.
template block2_gf< dlr_imfreq, matrix_valued > triqs::modest::gloc | ( | dlr_imfreq const & | mesh, |
one_body_elements_on_grid const & | obe, | ||
double | mu | ||
) |
template block2_gf< imfreq, matrix_valued > triqs::modest::gloc | ( | imfreq const & | mesh, |
one_body_elements_on_grid const & | obe, | ||
double | mu | ||
) |
template block2_gf< dlr_imfreq, matrix_valued > triqs::modest::gloc | ( | one_body_elements_on_grid const & | one_body, |
double | mu, | ||
block2_gf< dlr_imfreq, matrix_valued > const & | Sigma_dynamic, | ||
nda::array< nda::matrix< dcomplex >, 2 > const & | Sigma_static | ||
) |
template block2_gf< imfreq, matrix_valued > triqs::modest::gloc | ( | one_body_elements_on_grid const & | one_body, |
double | mu, | ||
block2_gf< imfreq, matrix_valued > const & | Sigma_dynamic, | ||
nda::array< nda::matrix< dcomplex >, 2 > const & | Sigma_static | ||
) |
void triqs::modest::h5_read | ( | h5::group | g, |
std::string const & | name, | ||
atomic_orbs & | shell | ||
) |
void triqs::modest::h5_read | ( | h5::group | g, |
std::string const & | name, | ||
band_dispersion & | bd | ||
) |
void triqs::modest::h5_read | ( | h5::group | g, |
std::string const & | name, | ||
downfolding_projector & | proj | ||
) |
void triqs::modest::h5_read | ( | h5::group | g, |
std::string const & | name, | ||
embedding & | embed | ||
) |
void triqs::modest::h5_read | ( | h5::group | g, |
std::string const & | name, | ||
initial_data & | data | ||
) |
Definition at line 112 of file checkpoint.hpp.
void triqs::modest::h5_read | ( | h5::group | g, |
std::string const & | name, | ||
iteration_data & | data | ||
) |
Definition at line 158 of file checkpoint.hpp.
void triqs::modest::h5_read | ( | h5::group | g, |
std::string const & | name, | ||
local_space & | ls | ||
) |
void triqs::modest::h5_read | ( | h5::group | g, |
std::string const & | name, | ||
one_body_elements_on_grid & | x | ||
) |
void triqs::modest::h5_read | ( | h5::group | g, |
std::string const & | name, | ||
spin_kind_e & | spin_kind | ||
) |
void triqs::modest::h5_write | ( | h5::group | g, |
std::string const & | name, | ||
atomic_orbs const & | shell | ||
) |
void triqs::modest::h5_write | ( | h5::group | g, |
std::string const & | name, | ||
band_dispersion const & | bd | ||
) |
void triqs::modest::h5_write | ( | h5::group | g, |
std::string const & | name, | ||
downfolding_projector const & | proj | ||
) |
void triqs::modest::h5_write | ( | h5::group | g, |
std::string const & | name, | ||
embedding const & | embed | ||
) |
void triqs::modest::h5_write | ( | h5::group | g, |
std::string const & | name, | ||
initial_data const & | data | ||
) |
Definition at line 105 of file checkpoint.hpp.
void triqs::modest::h5_write | ( | h5::group | g, |
std::string const & | name, | ||
iteration_data const & | data | ||
) |
Definition at line 146 of file checkpoint.hpp.
void triqs::modest::h5_write | ( | h5::group | g, |
std::string const & | name, | ||
local_space const & | ls | ||
) |
void triqs::modest::h5_write | ( | h5::group | g, |
std::string const & | name, | ||
one_body_elements_on_grid const & | x | ||
) |
void triqs::modest::h5_write | ( | h5::group | g, |
std::string const & | name, | ||
spin_kind_e const & | spin_kind | ||
) |
nda::array< nda::matrix< dcomplex >, 2 > triqs::modest::Hloc | ( | std::vector< tb_hamiltonian > const & | H_sigma, |
std::vector< atomic_orbs > const & | atomic_shells | ||
) |
Compute Hloc = H(R=0) given n_sigma tight_binding Hamiltonians.
Helper function to calculate Hloc given a vector of H with length sigma (contained in Wannier OBE)
H_sigma | a list of TB Hamiltonians of length n_sigma |
atomic_shells | a list of atomic shells corresponding to the orbitals contained in the TB Hamiltonian |
Definition at line 64 of file obe_tb.cpp.
template block_gf< dlr_imfreq, matrix_valued > triqs::modest::hybridization | ( | std::vector< nda::matrix< dcomplex > > const & | epsilon_levels, |
block_gf< dlr_imfreq, matrix_valued > const & | Gloc | ||
) |
template block_gf< dlr_imfreq, matrix_valued > triqs::modest::hybridization | ( | std::vector< nda::matrix< dcomplex > > const & | epsilon_levels, |
block_gf< dlr_imfreq, matrix_valued > const & | Gloc, | ||
block_gf< dlr_imfreq, matrix_valued > const & | Sigma_dynamic, | ||
std::vector< nda::matrix< dcomplex > > const & | Sigma_static | ||
) |
template block_gf< imfreq, matrix_valued > triqs::modest::hybridization | ( | std::vector< nda::matrix< dcomplex > > const & | epsilon_levels, |
block_gf< imfreq, matrix_valued > const & | Gloc | ||
) |
template block_gf< imfreq, matrix_valued > triqs::modest::hybridization | ( | std::vector< nda::matrix< dcomplex > > const & | epsilon_levels, |
block_gf< imfreq, matrix_valued > const & | Gloc, | ||
block_gf< imfreq, matrix_valued > const & | Sigma_dynamic, | ||
std::vector< nda::matrix< dcomplex > > const & | Sigma_static | ||
) |
triqs::modest::INSTANTIATE | ( | dlr_imfreq | ) |
triqs::modest::INSTANTIATE | ( | imfreq | ) |
triqs::modest::INSTANTIATE | ( | refreq | ) |
|
constexpr |
Definition at line 10 of file postprocess.cpp.
nda::array< dcomplex, 4 > triqs::modest::load_rotate_and_format_projectors | ( | std::string const & | filename, |
ReadMode | mode, | ||
std::vector< cmat_t > const & | rot_mats, | ||
std::vector< long > const & | atom_decomp | ||
) |
Read projectors using ReadMode, rotate to local frame, and embed them in the M space. (internal)
Definition at line 114 of file loaders.cpp.
embedding triqs::modest::make_embedding_impl | ( | local_space const & | C_space, |
nda::array< std::vector< long >, 2 > const & | block_decomposition, | ||
std::optional< std::vector< long > > const & | atom_to_imp | ||
) |
Definition at line 39 of file embedding.cpp.
|
inline |
Definition at line 379 of file embedding.hpp.
embedding triqs::modest::make_embedding_with_equivalences | ( | local_space const & | C_space | ) |
Definition at line 64 of file embedding.cpp.
embedding triqs::modest::make_embedding_with_no_equivalences | ( | local_space const & | C_space | ) |
Definition at line 80 of file embedding.cpp.
C2PY_IGNORE one_body_elements_tb triqs::modest::make_obe_from_tb | ( | std::vector< tb_hamiltonian > | H_sigma, |
spin_kind_e | spin_kind, | ||
std::vector< atomic_orbs > | atomic_shells | ||
) |
Helper to contruct and return an OBE_tb object given a list of tb_Hamiltonians of length n_sigma.
Definition at line 126 of file obe_tb.cpp.
one_body_elements_tb triqs::modest::one_body_elements_from_wannier90 | ( | std::string const & | wannier_file_path, |
spin_kind_e | spin_kind, | ||
std::vector< atomic_orbs > | atomic_shells | ||
) |
Construct a obe_tb from Wannier90 in the case of a single spin index.
wannier_file_path | string to Wannier90 files, including the prefix, as in "path/to/file/seedname" to specify a Wannier files named in the format "seedname_tb.dat" |
spin_kind | spin type for this calculation |
atomic_shells | list of atomic shells input by the user |
Definition at line 23 of file obe_tb.cpp.
one_body_elements_tb triqs::modest::one_body_elements_from_wannier90 | ( | std::string const & | wannier_file_path_up, |
std::string const & | wannier_file_path_dn, | ||
spin_kind_e | spin_kind, | ||
std::vector< atomic_orbs > | atomic_shells | ||
) |
Construct a obe_tb from Wannier90 in the case with separate spin up/spin down channels.
wannier_file_path_up | string to Wannier90 files, including the prefix, for the up spin channel, as in "path/to/file/seedname" to specify a Wannier files named in the format "seedname_tb.dat" |
wannier_file_path_dn | string to Wannier90 files, including the prefix, for the down spin channel as in "path/to/file/seedname" to specify a Wannier files named in the format "seedname_tb.dat" |
spin_kind | spin type for this calculation |
atomic_shells | list of atomic shells input by the user |
Definition at line 39 of file obe_tb.cpp.
std::ostream & triqs::modest::operator<< | ( | std::ostream & | out, |
band_dispersion const & | x | ||
) |
Definition at line 14 of file printing.cpp.
std::ostream & triqs::modest::operator<< | ( | std::ostream & | out, |
downfolding_projector const & | x | ||
) |
Definition at line 26 of file printing.cpp.
std::ostream & triqs::modest::operator<< | ( | std::ostream & | out, |
embedding const & | E | ||
) |
Definition at line 157 of file embedding.cpp.
std::ostream & triqs::modest::operator<< | ( | std::ostream & | out, |
ibz_symmetry_ops const & | ibz | ||
) |
Definition at line 63 of file printing.cpp.
std::ostream & triqs::modest::operator<< | ( | std::ostream & | out, |
local_space const & | bd | ||
) |
Definition at line 35 of file printing.cpp.
std::ostream & triqs::modest::operator<< | ( | std::ostream & | out, |
one_body_elements_on_grid const & | obe | ||
) |
Definition at line 73 of file printing.cpp.
one_body_elements_on_grid triqs::modest::permute_local_space | ( | std::vector< std::vector< long > > const & | atom_partition, |
one_body_elements_on_grid const & | obe | ||
) |
Definition at line 50 of file downfolding.cpp.
std::vector< atomic_orbs > triqs::modest::read_atomic_shells | ( | auto const & | filename, |
ReadMode | mode | ||
) |
Read atomic shells according to ReadMode. (internal)
Definition at line 179 of file loaders.cpp.
std::tuple< nda::array< dcomplex, 4 >, nda::matrix< long >, nda::array< double, 1 > > triqs::modest::read_bands_and_weights | ( | std::string | filename, |
ReadMode | mode | ||
) |
Read band dispersion and k-weights according to ReadMode. (internal)
Definition at line 155 of file loaders.cpp.
ibz_symmetry_ops triqs::modest::read_ibz_symmetry_ops | ( | auto const & | filename, |
ReadMode | mode | ||
) |
Construct the ibz_symmetry_ops according to ReadMode.
Definition at line 201 of file loaders.cpp.
std::vector< cmat_t > triqs::modest::read_rotation_matrices | ( | std::string const & | filename, |
ReadMode | mode | ||
) |
Read rotation matrices from hdf5. (internal)
Definition at line 91 of file loaders.cpp.
nda::array< nda::matrix< dcomplex >, 1 > triqs::modest::read_spherical_to_dft_basis | ( | std::string | dft, |
std::vector< atomic_orbs > const & | atomic_shells | ||
) |
Prepare the spherical Ylm to DFT orbital basis rotations. (internal)
Definition at line 191 of file loaders.cpp.
spin_kind_e triqs::modest::read_spin_kind | ( | auto const & | filename | ) |
Setup spin_kind enum. (internal)
Definition at line 171 of file loaders.cpp.
|
inline |
Change basis.
U | Rotations by block in atomic decomposition |
x | ibz_symmetry_ops to rotate |
Definition at line 96 of file ibz_symmetry_ops.hpp.
one_body_elements_on_grid triqs::modest::rotate_local_basis | ( | nda::array< nda::matrix< dcomplex >, 2 > const & | U, |
one_body_elements_on_grid const & | x | ||
) |
Rotates the local basis of the downfolding projector.
Definition at line 42 of file downfolding.cpp.
|
inline |
Definition at line 18 of file downfolding.hpp.
utility to flatten a nested vector (move to utils/ ?)
free function for svd (FIXME: nda::linalg::svd) (see TRIQS/nda PR#103)
Definition at line 17 of file loaders.cpp.