The Model Hamiltonian¶
- ggce.model.model_coupling_map(coupling_type, t, Omega, lam)[source]¶
Returns the value for g, the scalar that multiplies the coupling in the Hamiltonian. Converts the user-input lambda value to this g. Uses pre-defined values for the dimensionless coupling to get g for a variety of pre-defined default models.
- Parameters:
coupling_type (str) – The desired coupling type. Can be Holstein, Peierls, BondPeierls, or EdwardsFermionBoson.
t (float) – The hopping strength.
Omega (float) – The (Einstein) phonon frequency. Absolute value is taken, since it will be negative sometimes (in the case of TFD).
lam (float) – The dimensionless coupling.
- Returns:
The value for the coupling (g).
- Return type:
- Raises:
RuntimeError – If an unknown coupling type is provided.
- class ggce.model.SingleTerm(coupling_type, psi, phi, dag, coupling, phonon_index, phonon_frequency)[source]¶
Bases:
MSONable
Contains the shift indexes, dagger status, coupling strength and boson type of a single term. Particularly, this single term corresponds to a single element of the sum in equation (10) in this PRB,
\[g \sum_i c_i^\dagger c_{i + \psi} b^d_{i + \phi}\]- property psi¶
The shift term on \(c\).
- Return type:
numpy.array
- property phi¶
The shift term on \(b^d\).
- Return type:
numpy.array
- property dag¶
The “dagger status” of the term. “+” for creation operator and “-” for annihilation operator. Must be either “+” or “-“.
- Return type:
- property coupling¶
The coupling strength of the term. Directly multiplies the operators in the Hamiltonian. This is \(g\).
- Return type:
- class ggce.model.Hamiltonian(terms=[], phonon_frequencies=[], dimension=1, hopping=1.0)[source]¶
Bases:
MSONable
Container for all relevant parameters for defining the Hamiltonian. Current available are Holstein, EdwardsFermionBoson, BondPeierls, and Peierls.
- property terms¶
- property phonon_frequencies¶
- property dimension¶
- get_dict_rep()[source]¶
Provides a dictionary representation of the Hamiltonian. Packages terms by the phonon type.
- Return type:
- add_(coupling_type, phonon_index, phonon_frequency, coupling_strength=None, dimensionless_coupling_strength=None, coupling_multiplier=1.0)[source]¶
Adds a set of terms to the list of terms, corresponding to the provided coupling type, strength and index.
- Parameters:
coupling_type (str) – The type of coupling to add to the Hamiltonian.
phonon_index (int) – The index of the phonon to register.
phonon_frequency (float) – The Einstein phonon frequency corresponding to the phonon creation and annihilation operators. Traditionally denoted \(\Omega\).
coupling_strength (float) – The strength of the coupling to add. This is g. Default is None. If None, requires dimensionless_coupling_strength to be set.
dimensionless_coupling_strength (float) – The dimensionless coupling strength, which is a function of the hopping, and phonon frequency, as well as the type of model. Default is None. If None, requires coupling_strength to be set.
coupling_multiplier (float, optional) – Extra term which will multiply the coupling. This could be a prefactor from e.g. the TFD formalism. Default is 1.0
- class ggce.model.Model(lattice_constant=1.0, temperature=0.0, hamiltonian=, phonon_max_per_site=None, phonon_extent=[], phonon_number=[], phonon_absolute_extent=None, n_phonon_types=0)[source]¶
Bases:
MSONable
The Model object. Contains all information to be fed to the System class for solving for the Green’s function.
Warning
Note it is highly recommended to use the
from_parameters
classmethod instead, since the constructor is primarily designed to be used with MSONable and not all parameters are intended to be manually modified.- property temperature¶
The temperature for a TFD simulation, if requested (the default is 0.0).
- Return type:
- property phonon_absolute_extent¶
The absolute extent value. Controls how much the clouds can overlap. Indicates the largest possible cloud size when multi-phonon mode calculations are run.
- Return type:
- property n_phonon_types¶
A counter which keeps track of the number of phonon types that are contained in the model
- Return type:
- property phonon_max_per_site¶
Controls the hard-boson constraint. Essentially restrains the number of phonons that can be on each site. If None, there is no restriction.
- Return type:
- property hamiltonian¶
The Hamiltonian object for this model.
- Return type:
- classmethod from_parameters(hopping=1.0, lattice_constant=1.0, temperature=0.0, phonon_max_per_site=None, dimension=1)[source]¶
Convenience method for initializing the Model object from intuitive parameters. It is recommended that this classmethod is used over the constructor.
- Parameters:
hopping (float) – The hopping strength \(t\). Default is 1.0.
lattice_constant (float, optional) – Default is 1.0.
temperature (float, optional) – Note that setting this > 0 will engage the Thermofield Double method. Default is 0.0.
phonon_max_per_site (None, optional) – Used for hard phonon constraints. Default is None, indicating no limitation on the number of phonons per site.
dimension (int, optional) – The spatial dimension of the system to consider. Default is 1.
- add_(coupling_type, phonon_frequency, phonon_extent, phonon_number, phonon_extent_tfd=None, phonon_number_tfd=None, coupling_strength=None, dimensionless_coupling_strength=None, phonon_index_override=None)[source]¶
Adds an electron-phonon contribution to the model. Note that the user can override the boson index manually to construct more complicated models, such as single phonon-mode Hamiltonians but with multiple contributions to the coupling_strength term.
- Parameters:
coupling_type (str) – The type of coupling to add to the Hamiltonian.
phonon_frequency (float) – The Einstein phonon frequency corresponding to the phonon creation and annihilation operators. Traditionally denoted \(\Omega\).
phonon_extent (int) – Positive number for the max phonon extent for this term.
phonon_number (int) – Positive number for the max number of phonons for this term.
phonon_extent_tfd (int, optional) – Positive number for the max phonon extent for the TFD version of the term. Default is None.
phonon_number_tfd (int, optional) – Positive number for the max number of phonons for the TFD version of this term. Default is None.
coupling_strength (float, optional) – The strength of the coupling to add. This is g. Default is None. If None, requires dimensionless_coupling_strength to be set.
dimensionless_coupling_strength (float, optional) – The dimensionless coupling strength, which is a function of the hopping, and phonon frequency, as well as the type of model. Default is None. If None, requires coupling_strength to be set.
phonon_index_override (int, optional) – Overrides the internal phonon counter. Use at your own risk. Note that for finite-temperature simulations, these indexes should always come in pairs. In other words, if the override 0, then the TFD part of the terms will be 1. Throws an error if this number is odd and temperature > 0. Default is None (no override).
The Matrix Generation Engine¶
Terms¶
- ggce.engine.terms.config_legal(config, max_phonons_per_site=None, phonon_extent=None, phonon_number=None, allow_green=False)[source]¶
Helper method designed to test a standalone configuration of phonons for legality. A legal
config
array satisfies the following properties:It has at least two axes (one phonon index and one spatial)
It has only nonnegative entries
Every spatial “edge” has at least one phonon of some type
The maximum number of phonons per site condition is satisfied
The phonon extent criterion is satisfied for all phonon types
- Parameters:
config (numpy.array) – The input configuration.
max_phonons_per_site (int, optional) – If not None, checks that the config satisfies the maximum number of phonons per site. Returns False if any site contains more than the specified number.
phonon_extent (list, optional) – A list of int where each entry is the extent allowed for that phonon type. Checks this if not None, otherwise ignores.
allow_green (bool, optional) – If True, automatically passes the check if the total number of phonons sum to 0.
- Returns:
Returns True if the config passes all the tests. False if it fails any of them.
- Return type:
- class ggce.engine.terms.Config(config, max_modifications=1, modifications=0)[source]¶
Bases:
MSONable
A class for holding phonon occupations and defining operations on the cloud.
- property n_spatial_dimensions¶
The number of spatial dimensions in the config. Explicitly excludes counting the first dimension, which corresponds to the phonon index.
- Return type:
- property config¶
An array of the shape (n_phonon_types, cloud length axis 1, 2, …). The array should only contain integers. Spatial information is contained in the indexes 1 and greater. The 0th index contains the information about the phonon type.
- Return type:
- property n_phonon_types¶
The number of different types of phonons in the config. This is equal to the shape of the first axis of config.
- Return type:
- property phonon_cloud_shape¶
The shape corresponding to this particular cloud, ignoring information about phonon type.
- Return type:
- property total_phonons_per_type¶
Gets the total number of phonons of each type. Essentially, sums over all axes except the first.
- Return type:
- property total_phonons¶
Sums over the entire config. Returns the total number of phonons in the config.
- Return type:
- validate()[source]¶
Checks the config for validity and throws a
logger.critical
, which also kills the program if the validity checks fail. If any of the following fail, a critical is raised:If any element in the config has less than 0 phonons.
If any of the edges of the cloud has 0 phonons total.
Warning
The Green’s function, i.e. the case where there are actually 0 phonons on the cloud, will not raise any critical.
- id()[source]¶
Returns the string representation of the config. This is constructed by flattening the config, converting to a list and then converting that list to a string. In order to ensure any possible ambiguity issues, the shape of the initial config is also used as an identifier.
- Return type:
- remove_phonon_(*indexes)[source]¶
Executes the removal of a phonon, and then all following phonon removal “reduction rules”. I.e., removes slices of all zeros from the edges of the phonon clouds. Essentially Appendix A in this PRB. Specifically equations A1 and A3.
- Parameters:
*indexes – A tuple that must be of the same dimension as the config itself. The first value is the phonon index, and the others correspond to the real-space location in the phonon to remove. These indexes must fall within the existing config, else the program will abort.
- Returns:
An array of the index shifts used for modifying the coefficients of auxiliary Green’s functions.
- Return type:
- add_phonon_(*indexes)[source]¶
Executes the addition of a phonon, and then all following phonon addition “reduction rules”. Essentially Appendix A in this PRB. Specifically equations A2 and A4.
- Parameters:
*indexes – A tuple that must be of the same dimension as the config itself. The first value is the phonon index, and the others correspond to the real-space location in the phonon to remove. These indexes can fall outside of the current config, since phonons can be added anywhere on the chain.
- Returns:
An array of the index shifts used for modifying the coefficients of auxiliary Green’s functions.
- Return type:
- class ggce.engine.terms.Term(config, hamiltonian_term=None, constant_prefactor=1.0, exp_shift=None, f_arg=None, g_arg=None)[source]¶
Bases:
MSONable
Base class for a single term in an equation. This object contains all information required for defining a term in the sum of equation 22 in the PRB. that this code is based on.
- property config¶
The configuration of phonons. The Config class contains the necessary methods for adding and removing phonons of various types and in various locations. See
Config
.- Return type:
- property hamiltonian_term¶
A single term in the coupling of the Hamiltonian.
- Return type:
- property constant_prefactor¶
A constant prefactor that multiplies the evaluated value of the entire term when constructing the eventual matrix to invert at the end of the software pipeline.
- Return type:
- property exp_shift¶
Factor which multiplies the entire final prefactor term in addition to
Term.constant_prefactor
. It is given by\[e^{i \mathbf{k} \cdot \boldsymbol{\delta} \odot \mathbf{a}}\]where \(\delta\) is this
exp_shift
. Furthermore, this term is determined by the relations between the different f-functions. It is modified during the reduction rules as noted inConfig.remove_phonon_
andConfig.add_phonon_
.Warning
The shape of the
exp_shift
should be equal to the number of spatial dimensions of the phonon configuration.- Return type:
- property f_arg¶
This array corresponds to the value of the argument of the auxiliary Green’s function, \(f_n(x)\), where \(x\) is the
f_arg
. If None, this implies that thef_arg
value is left general. This should occur if a derived class isIndex
-like, meaning it indexes an entire equation. This is a term that originally appears on the “left hand side” of the equation of motion. During generation of the specific equations, this will be set to a “true delta” value and that equation appended to the master list.- Return type:
- property g_arg¶
Identical to
Term.f_arg
, but for the argument of the lattice Green’s function \(g_0(x)\), where \(x\) isg_arg
.- Return type:
- id(full=False)[source]¶
Returns a string with which one can index the term. There are two types of identifiers the user can request. First is the full=False version, which produces a string of the form
X(f_arg)
, whereX
is the configuration andf_arg
is the argument of the auxiliary Green’s function. Second is full=True, which returns the same string as full=False, but withg_arg
andc_exp
appended to the end as well.g_arg
is the argument of \(g_0\) andc_exp
is the exponential shift indexes.
- update_phonon_config_()[source]¶
Specific update scheme needs to be run depending on whether we add or remove a phonon.
- check_if_green_and_simplify_()[source]¶
There are cases when the algorithm produces results that look like e.g.
G(1)[0]
. This corresponds tof_0(1)
, which is actually just a Green’s function times a phase factor. The precise relation is\[f_0(\boldsymbol{\delta}) = e^{i \mathbf{k} \cdot \boldsymbol{\delta} \odot \mathbf{a}} G(k, w)\]We need to simplify this back to a term like G(0)[0], as if we don’t, when the specific equations are generated, it will think there are multiple Greens equations, of which of course there can be only one.
- class ggce.engine.terms.IndexTerm(config)[source]¶
Bases:
Term
A term that corresponds to an index of the equation. Critical differences between this and the base class include that the
f_arg
object is set to None by default and the constant_prefactor is set to 1.
- class ggce.engine.terms.EOMTerm(config, hamiltonian_term, model)[source]¶
Bases:
Term
A special instance of the base class that is part of the EOM of the Green’s function.
- coefficient(k, w, eta)[source]¶
This will set the prefactor to G0, since we assume that as a base class it originates from the g0 terms in the main EOM. Note that an index term does not have this method defined, since that terms prefactor should always be 1 (set by default to the constant prefactor), and this method will be overridden by AnnihilationTerm and CreationTerm classes.
- class ggce.engine.terms.NonIndexTerm(config, hamiltonian_term, model, constant_prefactor)[source]¶
Bases:
Term
- class ggce.engine.terms.AnnihilationTerm(config, hamiltonian_term, model, constant_prefactor)[source]¶
Bases:
NonIndexTerm
Specific object for containing f-terms in which we must subtract a phonon from a specified site.
- class ggce.engine.terms.CreationTerm(config, hamiltonian_term, model, constant_prefactor)[source]¶
Bases:
NonIndexTerm
Specific object for containing f-terms in which we must subtract a phonon from a specified site.
Equations¶
The equations
module constructs intermediate objects for dealing with
individual equations in the overall closure.
An Equation
can be loosely thought of in the following
pseudo-mathematical equation:
where Auxiliary Green’s functions, \(f\), are indexed by the total number of phonons, \(n\). Note that this is of course a “bad” quantum number and that much more information is required to represent the Auxiliary Green’s function, and that this information is actually contained, but for the sake of explanation we suppress those parameters.
For couplings of the forms compatible with the GGCE software package, only
a single phonon is created or annihilated at a time. Thus, the equations in
the closer couple Auxiliary Green’s Functions with \(n\) phonons to those
with \(n \pm 1\) phonons. There is also a bias term, \(b\), which is
generally 0 except in the case of the GreenEquation
, which is a
special instance of the Equation
used only for the true Green’s
function.
Tip
During the operation of the GGCE code, it is unlikely you will need to
actually manipulate the Equation
and GreenEquation
objects at all, so for the pure user, this module is likely not so
necessary to understand.
Note
See the Equation
and GreenEquation
for more details on
these two important objects.
- class ggce.engine.equations.Equation(index_term, model, f_arg_terms=None, terms_list=None)[source]¶
Bases:
MSONable
A single equation that indexes a single Auxiliary Green’s Function (AGF). These functions can be loosely conceptualized in terms of the number of phonons in each term, and can be written as
\[f(n) \sim \sum_i f_i(n-1) + \sum_j f_j(n+1) + b\]where \(n\) is an in-exhaustive index/bad quantum number corresponding to the total number of phonons on those sites. The
Equation
class contains anindex_term
, which in this case represents \(f(n)\), and aterms_list
object, which corresponds to the right hand side.There is also a bias term, \(b\), which represents possible inhomogeneities. Only the Green’s function itself contains \(b\neq 0\). All of these terms are implicitly functions of the momentum, \(k\), the energy/frequency, \(\omega\), and the artificial broadening \(\eta\). For the base
Equation
, thebias
is zero.Formally, this class is the representation of equation 22 in PRB. In one dimension, it is
\[f_\mathbf{n}(\delta) = \sum_{(g, \psi, \phi, \xi)} g \sum_{\gamma \in \Gamma_L^\xi} n^{(\xi, \gamma)} g_0(\delta + \gamma - \phi + \psi, \tilde{\omega}) f_{\mathbf{n}}^{(\xi, \gamma)}(\phi - \gamma)\](see the manuscript for a detailed explanation of all of these terms).
- classmethod from_config(config, model)[source]¶
Initializes the
Equation
class from a numpy array and aggce.models.Model
.- Parameters:
config (numpy.array) – The array representing the configuration of phonons.
model (ggce.models.Model) – Container for the full set of parameters. Contains all Hamiltonian and parameter information about the terms, phonons, couplings, etc.
- Return type:
- property index_term¶
The index term/config referencing the lhs of the equation.
- Return type:
- property model¶
Container for the full set of parameters. Contains all Hamiltonian and parameter information about the terms, phonons, couplings, etc.
- Returns:
Description
- Return type:
- property f_arg_terms¶
A dictionary of lists of ints indexed by the config string representation containing lists of integers highlighting the required values of delta for that particular config string representation.
- Return type:
- bias(k, w, eta)[source]¶
The default value for the bias is 0, except in the case of the Green’s function.
- visualize(full=True, coef=None)[source]¶
Prints the representative strings of the terms for both the indexer and terms. This allows the user to at a high level, and without coefficients, visualize all the terms in the equation. If coef is not None, we assume we want to visualize the coefficient too, and coef is a 2-vector containing the k and w points.
- Parameters:
full (bool, optional) – If True, prints out not only the configuration, its shape and the
f_arg
terms, but also theg_arg
andshift
components of theggce.engine.terms.Term
object. Default is True.coef (None, optional) – Optional coefficients for \(k\), \(\omega\) and \(\eta\) to pass to the terms and output their coefficients. This is primarily used for debugging. Default is None (indicating no coefficient information will be printed).
- class ggce.engine.equations.GreenEquation(model)[source]¶
Bases:
Equation
Equation object corresponding to the Green’s function. In one dimension, this corresponds directly with equation 13 in this PRB,
\[f_0(0) - \sum_{(g,\psi,\phi)} g e^{ik R_{\psi - \phi}} f_1(\phi) = G_0(k, \omega)\]The primary way in which this object differs from
Equation
is that it has a non-zeroGreenEquation.bias
.
System¶
- ggce.engine.system.config_space_gen(length, total_sum)[source]¶
Generator for yielding all possible combinations of integers of length
length
that sum tototal_sum
.Warning
Not that cases such as
length == 4
andtotal_sum == 5
such as[0, 0, 2, 3]
still need to be screened out, since these do not correspond to valid f-functions.Note
The algorithm to produce this code can be found here.
- ggce.engine.system.generate_all_legal_configurations(model)[source]¶
Summary
In one dimension, this is really easy. We can simply iterate over all 0 < n <= N and 0 < m <= M. Things get much more complicated as the dimensionality increases. However, it is simplified somewhat since we always assume that M is a “square”. We can still generate configurations that
- Parameters:
model (TYPE) – Description
- class ggce.engine.system.System(model=None, generalized_equations=None, f_arg_list=None, equations=None, root=None, mpi_comm=None, chkpt_lim=1000000.0, autoprime=True)[source]¶
Bases:
object
Defines a list of Equations (a system of equations, so to speak) and operations on that system.
- property model¶
- property generalized_equations¶
A list of equations constituting the system; in general form, meaning all except for the Green’s function are not defined for specific delta values.
- Returns:
A list of
Equation
objects.- Return type:
- property equations¶
- checkpoint()[source]¶
Runs the checkpoint protocol to attempt and save the current state of the system. This method will only run if a root directory is defined at class instantiation. Note that existing checkpoints will never be overwritten.
If the number of equations is larger than a pre-specified cutoff (by default 1,000,000), then checkpointing raises a warning and does not work.
- classmethod from_checkpoint(root)[source]¶
Reloads the state of the System from the saved pickle file checkpoint.
- Parameters:
root (os.PathLike) – Checkpoint directory
- visualize(generalized=True, full=True, coef=None)[source]¶
Allows for easy visualization of the closure. Note this isn’t recommended when there are greater than 10 or so equations, since it will be very difficult to see everything.
- Parameters:
generalized (bool) – If True, prints information on the generalized equations, else prints the full equations. Default is True.
full (bool) – If True, prints information about the argument of g(…) and the exponential shift in addition to the n-vector and f-argument. Else, just prints the latter two. Default is True.
coef (list, optional) – If not None, actually evaluates the terms at the value of the (k, w) coefficient. Default is None.
- get_basis(full_basis=False)[source]¶
Prepares the solver-specific information.
Returns the non-zero elements of the matrix in the following format. The returned quantity is a dictionary indexed by the order of the hierarchy (in this case, the number of phonons contained). Each element of this dictionary is another dictionary, with the keys being the index term identifier (basically indexing the row of the matrix), and the values a list of tuples, where the first element of each is the identifier (a string) and the second is a callable function of k and omega, representing the coefficient at that point.
Solvers¶
- class ggce.executors.solvers.Solver(system=None, root=None, basis=None, mpi_comm=None)[source]¶
Bases:
ABC
- property system¶
- property root¶
- property basis¶
- property mpi_comm¶
- property mpi_rank¶
- property mpi_world_size¶
- class ggce.executors.solvers.BasicSolver(system=None, root=None, basis=None, mpi_comm=None)[source]¶
Bases:
Solver
- class ggce.executors.solvers.SparseSolver(*args, **kwargs)[source]¶
Bases:
BasicSolver
A sparse, serial solver for the Green’s function. Useful for when the calculation being performed is quite cheap. Note that there are a variety of checkpointing features automatically executed when paths are provided.
When a
System
and path are provided, the system’s root directory is overwritten and that object immediately checkpointed to the provided directory.If only a
System
object is provided, no checkpointing will be performed.If only a path is provided, the solver will attempt to instantiate the
System
object. Checkpointing will the proceed as normal.
- class ggce.executors.solvers.DenseSolver(*args, **kwargs)[source]¶
Bases:
BasicSolver
A sparse, dense solver for the Green’s function. Note that there are a variety of checkpointing features automatically executed when paths are provided. Uses the SciPy dense solver engine to solve for G(k, w) in serial. This method uses the continued fraction approach,
\[R_{n-1} = (1 - \beta_{n-1}R_{n})^{-1} \alpha_{n-1}\]with
\[R_n = \alpha_n\]and
\[f_n = R_n f_{n-1}\]