Parquet Algorithm to Build Diagrams

API

FeynmanDiagram.FrontEnds.Parquet.ParquetBlocksType
struct ParquetBlocks

The channels of the left and right sub-vertex4 of a bubble diagram in the parquet equation

Members

  • phi : channels of left sub-vertex for the particle-hole and particle-hole-exchange bubbles
  • ppi : channels of left sub-vertex for the particle-particle bubble
  • Γ4 : channels of right sub-vertex of all channels
source
FeynmanDiagram.FrontEnds.Parquet.ep_couplingMethod
function ep_coupling(para::DiagPara;
    extK=[getK(para.totalLoopNum, 1), getK(para.totalLoopNum, 2), getK(para.totalLoopNum, 3)],
    channels::AbstractVector=[PHr, PHEr, PPr, Alli],
    subdiagram=false,
    name=:none, resetuid=false,
    blocks::ParquetBlocks=ParquetBlocks()
)

Generate electron-phonon 4-vertex diagrams using Parquet Algorithm. The right incoming Tau will be set to the last Tau for all diagrams | | Γ3 –––-| | |

Arguments

  • para : parameters. It should provide internalLoopNum, interactionTauNum, firstTauIdx
  • extK : basis of external loops as a vector [left in, left out, right in, right out].
  • channels : vector of channels in the left Γ3 diagrams.
  • subdiagram : a sub-vertex or not
  • name : name of the vertex
  • resetuid : restart uid count from 1
  • blocks : building blocks of the Parquet equation. See the struct ParquetBlocks for more details.

Output

  • A DataFrame with fields :response, :extT, :diagram, :hash.

Output

  • A DataFrame with fields :response, :type, :extT, :diagram, :hash
source
FeynmanDiagram.FrontEnds.Parquet.get_ver4IMethod
get_ver4I()

Retrieves the global dictionary `vertex4I_diags` that contains graph initialized by `initialize_vertex4I_diags`. 
This function is a getter that provides access to the stored graph data of the 3- and 4-point fully-irreducible (Alli) vertex functions.
source
FeynmanDiagram.FrontEnds.Parquet.greenFunction
green(para::DiagPara, extK = getK(para.totalLoopNum, 1), extT = para.hasTau ? (1, 2) : (0, 0), subdiagram = false;
    name = :G, resetuid = false, blocks::ParquetBlocks=ParquetBlocks())

Build composite Green's function. By definition, para.firstTauIdx is the first Tau index of the left most self-energy subdiagram.

Arguments

  • para : parameters. It should provide internalLoopNum, interactionTauNum, firstTauIdx
  • extK : basis of external loop.
  • extT: [Tau index of the left leg, Tau index of the right leg]
  • subdiagram : a sub-vertex or not
  • name : name of the diagram
  • resetuid : restart uid count from 1
  • blocks : building blocks of the Parquet equation. See the struct ParquetBlocks for more details.

Output

  • A Graph object or nothing if the Green's function is illegal.
source
FeynmanDiagram.FrontEnds.Parquet.initialize_vertex4I_diagsMethod
initialize_vertex4I_diags(; filter=[NoHartree], spinPolarPara::Float64=0.0)

Initialize the vertex4I_diags dictionary with the diagrams of the 3- and 4-point fully-irreducible (Alli) vertex functions.

Parameters

  • filter (optional) : a list of filter conditions to select the diagrams. Default is [NoHartree].
  • spinPolarPara (optional) : the spin-polarization parameter. Default is 0.0.
source
FeynmanDiagram.FrontEnds.Parquet.innerTauNumMethod
function innerTauNum(type::DiagramType, innerLoopNum, interactionTauNum)

internal imaginary-time degrees of freedom for a given diagram type and internal loop number.
For the vertex functions (self-energy, polarization, vertex3, and vertex4), innerTauNum is equivalent to tauNum.
For the Green function, tauNum = innerTauNum + external tauNum
source
FeynmanDiagram.FrontEnds.Parquet.mergebyFunction
function mergeby(df::DataFrame, fields=Vector{Symbol}();
    operator=Sum(), name::Symbol=:none,
    getid::Function=g -> GenericId(g[1, :diagram].properties.para, Tuple(g[1, fields]))
)

Aggregates and merges rows in a DataFrame based on specified fields and an aggregation operator. It is designed to work with data frames containing diagram representations or similar structured data, allowing for flexible aggregation based on custom identifiers and properties.

Parameters

  • df: A DataFrame to be processed. It must contain a column named :diagram which holds Graph.
  • fields: A vector of Symbols specifying the columns based on which the aggregation groups are formed.
  • operator: The aggregation operator to be applied. This operator is used to aggregate data across rows within each group formed by the specified fields. The default is Sum().
  • name: An optional Symbol representing the name of a new column to store the aggregated results. The default is :none.
  • getid: A function that generates a unique identifier for each group based on the specified fields and the properties of the :diagram column. The default function is g -> GenericId(g[1, :diagram].properties.para, Tuple(g[1, fields])).
source
FeynmanDiagram.FrontEnds.Parquet.polarizationFunction
function polarization(para::DiagPara, extK=getK(para.totalLoopNum, 1), subdiagram=false; name=:Π, resetuid=false, blocks::ParquetBlocks=ParquetBlocks())

Generate polarization diagrams using Parquet Algorithm.

Arguments

  • para : parameters. It should provide internalLoopNum, interactionTauNum, firstTauIdx
  • extK : basis of external loop.
  • subdiagram : a sub-vertex or not
  • name : name of the vertex
  • resetuid : restart uid count from 1
  • blocks : building blocks of the Parquet equation. See the struct ParquetBlocks for more details.

Output

  • A DataFrame with fields :response, :diagram, :hash.
  • All polarization share the same external Tau index. With imaginary-time variables, they are extT = (para.firstTauIdx, para.firstTauIdx+1)
source
FeynmanDiagram.FrontEnds.Parquet.sigmaFunction
function sigma(para::DiagPara, extK=getK(para.totalLoopNum, 1), subdiagram=false;
    name=:Σ, resetuid=false, blocks::ParquetBlocks=ParquetBlocks())

Build sigma diagram. When sigma is created as a subdiagram, then no Fock diagram is generated if para.filter contains NoFock, and no sigma diagram is generated if para.filter contains Girreducible

Arguments

  • para : parameters. It should provide internalLoopNum, interactionTauNum, firstTauIdx
  • extK : basis of external loop.
  • subdiagram : a sub-vertex or not
  • name : name of the diagram
  • resetuid : restart uid count from 1
  • blocks : building blocks of the Parquet equation. See the struct ParquetBlocks for more details.

Output

  • A DataFrame with fields :type, :extT, :diagram, :hash
  • All sigma share the same incoming Tau index, but not the outgoing one
source
FeynmanDiagram.FrontEnds.Parquet.update_extKTFunction
update_extKT(diags::Vector{Graph}, para::DiagPara, legK::Vector{Vector{Float64}}) -> Vector{Graph}

Returns a new vector of graphs with updated external momenta (extK) and external times (extT), based on the provided graphs, parameters, and external legs' momenta.

Arguments

  • diags::Vector{Graph}: A vector of Graph objects.
  • para::DiagPara: parameters reconstructed in the graphs. Its firstTauIdx will update the extT of graphs.
  • legK::Vector{Vector{Float64}}: basis of the external momenta for the legs of the diagram as [left in, left out, right in, right out].
  • extraLoopIdx: the index of the extra loop in the external momenta basis in legK. Defaults to nothing, which means no extra loop is included.

Returns

  • Vector{Graph}: A new vector of Graph objects with updated extK, extT, and para (if existing) properties for each node.
source
FeynmanDiagram.FrontEnds.Parquet.update_extKT!Function
update_extKT!(diags::Vector{Graph}, para::DiagPara, legK::Vector{Vector{Float64}})

Update the external momenta (extK) and external times (extT) of all the nodes in a vector of graphs in-place.

Arguments

  • diags::Vector{Graph}: A vector of Graph objects.
  • para::DiagPara: parameters reconstructed in the graphs. Its firstTauIdx will update the extT of graphs.
  • legK::Vector{Vector{Float64}}: basis of the external momenta for the legs of the diagram as [left in, left out, right in, right out].
  • extraLoopIdx: the index of the extra loop in the external momenta basis in legK. Defaults to nothing, which means no extra loop is included.
source
FeynmanDiagram.FrontEnds.Parquet.vertex3Function
function vertex3(para::DiagPara, _extK=[getK(para.totalLoopNum, 1), getK(para.totalLoopNum, 2)], subdiagram=false;
    name=:Γ3, channels=[PHr, PHEr, PPr, Alli], resetuid=false, blocks::ParquetBlocks=ParquetBlocks())

Generate 3-vertex diagrams using Parquet Algorithm. With imaginary-time variables, all vertex3 generated has the same bosonic Tidx $extT[1]=para.firstTauIdx$ and the incoming fermionic Tidx $extT[2]=para.firstTauIdx+1$.

#Arguments

  • para : parameters. It should provide internalLoopNum, interactionTauNum, firstTauIdx
  • extK : basis of external loops as a vector [bosonic leg (out), fermionic in, fermionic out], extK[1] = extK[2] - extK[3].
  • subdiagram : a sub-vertex or not
  • name : name of the vertex
  • channels : vector of channels of the current 4-vertex.
  • resetuid : restart uid count from 1
  • blocks : building blocks of the Parquet equation. See the struct ParquetBlocks for more details.

Output

  • A DataFrame with fields :response, :extT, :diagram, :hash.
source
FeynmanDiagram.FrontEnds.Parquet.vertex4Function
vertex4(para::DiagPara,
    extK = [getK(para.totalLoopNum, 1), getK(para.totalLoopNum, 2), getK(para.totalLoopNum, 3)],
    subdiagram = false;
    channels::AbstractVector = [PHr, PHEr, PPr, Alli],
    level = 1, name = :none, resetuid = false,
    blocks::ParquetBlocks=ParquetBlocks(),
    blockstoplevel::ParquetBlocks=blocks
    )

Generate 4-vertex diagrams using Parquet Algorithm

Arguments

  • para : parameters. It should provide internalLoopNum, interactionTauNum, firstTauIdx
  • extK : basis of external loops as a vector [left in, left out, right in, right out].
  • subdiagram : a sub-vertex or not
  • channels : vector of channels of the current 4-vertex.
  • name : name of the vertex
  • level : level in the diagram tree
  • resetuid : restart uid count from 1
  • blocks : building blocks of the Parquet equation. See the struct ParquetBlocks for more details.
  • blockstoplevel : building blocks of the Parquet equation at the toplevel. See the struct ParquetBlocks for more details.

Output

  • A DataFrame with fields :response, :type, :extT, :diagram, :hash
source

Usage

julia> using FeynmanDiagram
julia> para = DiagPara(type = Ver4Diag, innerLoopNum = 1, hasTau = true);ERROR: UndefVarError: `Ver4Diag` not defined in `Main` Suggestion: check for spelling errors or missing imports.
julia> Parquet.vertex4(para)ERROR: UndefVarError: `para` not defined in `Main` Suggestion: check for spelling errors or missing imports.
julia> para = DiagPara(type = Ver3Diag, innerLoopNum = 1, hasTau = true);ERROR: UndefVarError: `Ver3Diag` not defined in `Main` Suggestion: check for spelling errors or missing imports.
julia> Parquet.vertex3(para)ERROR: UndefVarError: `para` not defined in `Main` Suggestion: check for spelling errors or missing imports.
julia> para = DiagPara(type = SigmaDiag, innerLoopNum = 1, hasTau = true);ERROR: UndefVarError: `SigmaDiag` not defined in `Main` Suggestion: check for spelling errors or missing imports.
julia> Parquet.sigma(para)ERROR: UndefVarError: `para` not defined in `Main` Suggestion: check for spelling errors or missing imports.
julia> para = DiagPara(type = PolarDiag, innerLoopNum = 1, hasTau = true);ERROR: UndefVarError: `PolarDiag` not defined in `Main` Suggestion: check for spelling errors or missing imports.
julia> Parquet.polarization(para)ERROR: UndefVarError: `para` not defined in `Main` Suggestion: check for spelling errors or missing imports.