Struct lattice_model

Nested Relationships

Nested Types

Struct Documentation

struct lattice_model

description of the lattice model

Public Functions

lattice_model()

default constructor

double Potthoff_functional()
int neighbor_index(vector3D<int64_t> &R)

Finds the index of a neighboring super unit cell, given a position R.

matrix<Complex> periodize(const vector3D<double> &k, matrix<Complex> &mat)

Periodizes a matrix (Green function, CPT Green function or self-energy) into the spin-orbital matrix gn (of dimension nband*n_mixed)

Parameters:
  • k – [in] wavevector in the superdual basis

  • mat – [in] input matrix in full orbital space (Green function, CPT Green function or self-energy)

Returns:

a complex-valued, periodized matrix

pair<string, int> name_and_label(string &S, bool cut = false)

From a string, extract the name of the operator and the system label (the latter separated from the name by ‘_’)

Parameters:
  • S – [in] input string

  • cut – [in] if true, cuts the cluster suffix of the name

Returns:

a pair <name, label>

size_t to_GF_index(size_t I, size_t spin, bool anomalous, bool spin_down)

converts an index as used in lattice_operator::elements to an Green function index

Parameters:
  • I – index in lattice_operator

  • anomalous – true if the index pertains to a nambu index

  • spin_down – true if the index pertains to a spin down degree of freedom (in case of HS_mixing::up_down)

Returns:

the index

template<typename T>
void build_cluster_operators(lattice_operator &op)

Builds the cluster operators deriving from the lattice operators.

Parameters:

op – [in] lattice operator

vector<Complex> periodize(const vector3D<double> &k, vector<Complex> &mat)

Periodizes a vector into the spin-band vector gn (output, of dimension nband*n_mixed)

Parameters:
  • k – [in] wavevector in the dual supersubtype basis (superdual)

  • mat – [in] input vector in full orbital space

Returns:

a complex-valued, periodized vector

void add_anomalous_elements(vector<lattice_matrix_element> &E, int s1, int s2, int ni, int ni_opp, Complex z, latt_op_type SC)

add anomalous elements associated with a given pair of sites

Parameters:
  • E[in] : vector of matrix elements

  • s1[in] : site no 1

  • s2[in] : site no 2

  • ni[in] : neighbor index

  • ni_opp[in] : opposite neighbor index

  • z[in] : amplitude of the element

  • SC[in] : type of pairing

void add_chemical_potential()

Adds the chemical potential to the model (this operator is always present)

void anomalous_operator(const string &name, vector3D<int64_t> &link, complex<double> amplitude, int orb1, int orb2, const string &type)

Defines a lattice anomalous operator.

Parameters:
  • name – [in] name of the operator

  • link – [in] bond vector

  • amplitude – [in] overall multiplier

  • b1 – [in] lattice orbital label of the first electron

  • b2 – [in] lattice orbital label of the second electron

  • SC_str – [in] type of superconductor (singlet, dz, dx, dy)

void close_model(bool force = false)

consolidates some model data after reading the operators

void density_wave(const string &name, vector3D<int64_t> &link, complex<double> amplitude, int orb, vector3D<double> Q, double phase, const string &type)

Defines a density wave operator.

Parameters:
  • name – [in] name of the operator

  • cdw_link – [in] bond for the density wave

  • amplitude – [in] overall multiplier of the operator

  • orb – [in] lattice orbital label affected

  • Q – [in] wavevector

  • phase – [in] phase as it appears in the complex exponential

  • type – [in] type of density wave

void explicit_operator(const string &name, const string &type, const vector<tuple<vector3D<int64_t>, vector3D<int64_t>, complex<double>>> &elem, int tau, int sigma)

Defines a one-body operator explicitly.

Parameters:
  • name – [in] name of the operator

  • type – [in] type or operator (see enumerate latt_op_type)

  • elem – [in] list of elements specified by position label, link and amplitude.

  • tau – [in] same meaning as for hopping terms.

  • sigma – [in] same meaning as for hopping terms.

matrix<Complex> compact_tiling(const matrix<Complex> &A, const vector3D<double> &k)

Compact-tiling of a dim_GF x dim_GF matrix at wavevector k.

For each displacement d, all intra-cluster elements A[s,s’] with the same d are averaged (step 1), then inter-cluster elements A[s, f(s+d)] are added with the Bloch phase exp(i*k*R*2*pi) for the wrapping vector R (step 2). Uses tiling_data pre-computed by prepare_tiling_data().

Parameters:
  • A – [in] matrix to compact-tile, of size dim_GF x dim_GF

  • k – [in] wavevector in the superdual basis (same convention as M.k and periodize())

Returns:

the compact-tiled matrix Act of size dim_GF x dim_GF

std::optional<neighbor_match> find_second_site(int s1, const vector3D<int64_t> &link)

finds the second site, given a first site and a link

Parameters:
  • s1 – [in] site

  • link – [in] link

Returns:

the matching site (s2) and neighbor indices, or std::nullopt if there is no site at that position (e.g. graphene)

void neighbor_census()

Adds all possible neighbor vectors to the neighbor list by examining every position difference r2-r1 between sites of the super unit cell and applying that displacement to every site in the super unit cell.

Algorithm: (1) Loop over sites i1 (position r1) of the super unit cell. (2) Loop over sites i2 (position r2) of the super unit cell. (3) For each site s of the super unit cell, call find_second_site with link = r2-r1: this folds the displaced position back into the SUC, computes the neighbor vector, and registers it via neighbor_index if new.

void prepare_tiling_data()

Builds tiling_data: for each distinct spatial displacement d between intra-cluster site pairs, stores the list of (s1, s2, n) triplets where s1, s2 are site indices and n is the neighbor label (n=0 for intra-cluster, n>0 for inter-cluster wrapping).

The intra entries (n=0) are placed first; n_intra records their count.

void hopping_operator(const string &name, vector3D<int64_t> &link, double amplitude, int orb1, int orb2, int tau, int sigma)

Defines a lattice hopping operator.

Parameters:
  • name – [in] name of the operator

  • link – [in] bond vector

  • amplitude – [in] overall multiplier

  • b1 – [in] lattice orbital label of the first electron

  • b2 – [in] lattice orbital label of the second electron

  • tau – [in] Pauli matrix for the orbital part

  • sigma – [in] Pauli matrix for the spin part

void current_operator(const string &name, vector3D<int64_t> &link, double amplitude, int b1, int b2, int dir, bool pau = true)

Defines a lattice current operator.

Parameters:
  • name – [in] name of the operator

  • link – [in] bond vector

  • amplitude – [in] overall multiplier

  • b1 – [in] lattice orbital label of the first electron

  • b2 – [in] lattice orbital label of the second electron

  • dir – [in] direction of the current (0,1,2) for (x,y,z)

  • pau – [in] = true for current from real hopping, false for current from imaginary hopping

void interaction_operator(const string &name, vector3D<int64_t> &link, double amplitude, int orb1, int orb2, const string &type)

Defines a lattice interaction operator.

Parameters:
  • name – [in] name of the operator

  • link – [in] bond vector

  • amplitude – [in] overall multiplier

  • b1 – [in] lattice orbital label of the first electron

  • b2 – [in] lattice orbital label of the second electron

  • type – [in] type of interaction (Hubbard, i.e., density-density, if nothing is specified)

void one_body_matrix(lattice_operator &op)

builds the GF matrix elements of the lattice operator op

Parameters:

the[in] lattice operator

void post_parameter_consolidate(size_t label)

consolidates some model data after reading the first set of parameters

Parameters:

label – [in] label of the model instance used

void pre_operator_consolidate()

Finalizes some initialization of the lattice and superlattice, positions, etc.

void print(ostream &fout, bool asy_operators = false, bool asy_labels = false, bool asy_orb = false, bool asy_neighbors = false, bool asy_working_basis = false)

prints info about the model in a file (for debugging)

Parameters:
  • fout – [in] output channel

  • asy_operators – [in] if true, prints an asymptote file for each operator

  • bool – asy_labels [in] if true, adds labels to the asymptote graphics

  • bool – asy_orb [in] if true, adds lattice orbital labels to the asymptote graphics

  • bool – asy_neighbors [in] if true, adds the neighboring clusters to the asymptote graphics

  • bool – asy_working_basis [in] if true, plots in the working basis instead of the physical basis

matrix<Complex> lattice_hybridization(int iw, int ik)

extracts a particular hybridization matrix

extracts a particular lattice hybridization matrix

Public Members

basis3D dual

dual of lattice

basis3D phys

basis of the physical lattice

basis3D physdual

dual of basis of the physical lattice

basis3D superdual

dual of superlattice

bool bath_exists

true if at least one cluster has a bath of uncorrelated orbitals (DMFT)

bool is_closed

true if operators can no longer be added, as the first instance of the model was created

int mixing

spin-Nambu mixing of the model (cannot be changed once set)

int n_mixed

= 2 for mixing=1, 2 and 4 for mixing=3

int spatial_dimension

dimension of the superlattice

int nsys

total number of systems in the model

lattice3D superlattice

superlattice (for the super unit cell)

lattice3D unit_cell

lattice (useful in case of multiband models)

map<string, shared_ptr<lattice_operator>> term

pointers to the terms of the Hamiltonian

map<vector3D<int64_t>, size_t> folded_position_map

maps folded positions to site indices

map<vector3D<int64_t>, size_t> position_map

maps positions to site indices

size_t dim_GF

dimension of the Green function

size_t dim_reduced_GF

dimension of the reduced (periodized) Green function

size_t Lc

number of lattice sites in the repeated unit = sites.size()/n_band

size_t n_band

number of bands in the model

string name

name of the model

vector<cluster> clusters

list of clusters

vector<System> systems

list of systems

vector<int> inequiv

labels of the non-equivalent clusters

vector<int> is_nambu_index
vector<pair<int, int>> bonds
vector<shared_ptr<lattice_operator>> one_body_ops

pointers to the one-body or anomalous operators (for computing averages)

vector<site> sites

list of sites

vector<size_t> s_
vector<size_t> GF_dims
vector<size_t> GF_offset
vector<size_t> opposite_neighbor

indices of opposite neighbors

vector<size_t> reduced_Green_index
vector<vector<vector<Complex>>> pauli

Pauli matrices.

vector<vector3D<double>> Green_to_position
vector<vector3D<int64_t>> neighbor

list of neighboring superclusters

vector<tile_group> tiling_data

tiling groups, built once by prepare_tiling_data()

shared_ptr<parameter_set> param_set

parameter set structure

vector<string> sector_strings

Hilbert space sectors to explore in the different systems.

string hybrid_file

name of HDF5 file containing the external hybridization

shared_ptr<lattice_hybrid> hybrid

external hybridization

Public Static Attributes

static bool model_consolidated = false
struct tile_group

pre-computed connectivity for compact_tiling(): one group per distinct spatial displacement

Public Members

size_t n_intra

number of intra-cluster (n=0) entries at start of entries

vector<array<size_t, 3>> entries

{s1, s2, n}: spatial site indices and neighbor label