BMBPT Diagram

Routines and class for Bogoliubov MBPT diagrams.

class adg.bmbpt.BmbptFeynmanDiagram(nx_graph, tag_num)[source]

Bases: adg.diag.Diagram

Describes a BMBPT Feynman diagram with its related properties.

two_or_three_body

The 2 or 3-body characted of the vertices.

Type:int
time_tag

The tag number associated to the diagram’s associated TSD.

Type:int
tsd_is_tree

The tree or non-tree character of the associated TSD.

Type:bool
feynman_exp

The Feynman expression associated to the diagram.

Type:str
diag_exp

The Goldstone expression associated to the diagram.

Type:str
vert_exp

The expression associated to the vertices.

Type:list
hf_type

The Hartree-Fock, non-Hartree-Fock or Hartree-Fock for the energy operator only character of the graph.

Type:str
attribute_expressions(time_diag)[source]

Attribute the correct Feynman and Goldstone expressions.

Parameters:time_diag (TimeStructureDiagram) – The associated TSD.
attribute_qp_labels()[source]

Attribute the appropriate qp labels to the graph’s propagators.

extract_integral()[source]

Return the integral part of the Feynman expression of the diag.

Returns:The integral part of its Feynman expression.
Return type:(str)
extract_numerator()[source]

Return the numerator associated to a BMBPT graph.

Returns:The numerator of the graph.
Return type:(str)
has_crossing_sign()[source]

Return True for a minus sign associated with crossing propagators.

Use the fact that all lines propagate upwards and the canonical representation of the diagrams and vertices.

Returns:
Encode for the sign factor associated with crossing
propagators.
Return type:(bool)
multiplicity_symmetry_factor()[source]

Return the symmetry factor associated with propagators multiplicity.

Returns:The symmetry factor associated with equivalent lines.
Return type:(str)
time_tree_denominator(time_graph)[source]

Return the denominator for a time-tree graph.

Parameters:time_graph (NetworkX MultiDiGraph) – Its associated time-structure graph.
Returns:The denominator of the graph.
Return type:(str)
vertex_exchange_sym_factor()[source]

Return the symmetry factor associated with vertex exchange.

Returns:The symmetry factor for vertex exchange.
Return type:(str)
vertex_expression(vertex)[source]

Return the expression associated to a given vertex.

Parameters:vertex (int) – The vertex of interest in the graph.
write_diag_exps(latex_file, norder)[source]

Write the expressions associated to a diagram in the LaTeX file.

Parameters:
  • latex_file (file) – The LaTeX outputfile of the program.
  • norder (int) – The order in BMBPT formalism.
write_graph(latex_file, directory, write_time)[source]

Write the BMBPT graph and its associated TSD to the LaTeX file.

Parameters:
  • latex_file (file) – The LaTeX output file of the program.
  • directory (str) – The path to the result folder.
  • write_time (bool) – True if we want informations on the associated TSDs.
write_section(result, commands, diags_nbs)[source]

Write section and subsections for BMBPT result file.

Parameters:
  • result (file) – The LaTeX output file of the program.
  • commands (dict) – The flags associated with run management.
  • diags_nbs (dict) – The number of diagrams per type.
write_tsd_info(diagrams_time, latex_file)[source]

Write info related to the BMBPT associated TSD to the LaTeX file.

Parameters:
  • diagrams_time (list) – The associated TSDs.
  • latex_file (file) – The LaTeX output file of the program.
write_vertices_values(latex_file, mapping)[source]

Write the qp energies associated to each vertex of the diag.

Parameters:
  • latex_file (file) – The LaTeX output file of the program.
  • mapping (dict) – A mapping between the vertices in the diagram and the vertices in its euivalent TSD, since permutations between vertices are possible.
adg.bmbpt.check_unconnected_spawn(matrices, max_filled_vertex, length_mat)[source]

Exclude some matrices that would spawn unconnected diagrams.

Parameters:
  • matrices (list) – The adjacency matrices to be checked.
  • max_filled_vertex (int) – The furthest vertex until which the matrices have been filled.
  • length_mat (int) – The size of the square matrices.
>>> mats = [[[0, 2, 0], [2, 0, 0], [0, 0, 0]],                 [[0, 2, 1], [2, 0, 1], [0, 0, 0]]]
>>>
>>> check_unconnected_spawn(mats, 1, 3)
>>> mats
[[[0, 2, 1], [2, 0, 1], [0, 0, 0]]]
adg.bmbpt.diagrams_generation(p_order, three_body_use, nbody_obs, canonical)[source]

Generate diagrams for BMBPT from bottom up.

Parameters:
  • p_order (int) – The BMBPT perturbative order of the studied diagrams.
  • three_body_use (bool) – Flag for the use of three-body forces.
  • nbody_obs (int) – N-body character of the obervable of interest.
  • canonical (bool) – True if one draws only canonical diagrams.
Returns:

NumPy arrays encoding the adjacency matrices of the graphs.

Return type:

(list)

>>> diagrams_generation(1, False, 2, False) #doctest: +NORMALIZE_WHITESPACE
[array([[0, 4], [0, 0]]), array([[0, 2], [0, 0]])]
>>> diagrams_generation(1, True, 3, False)  #doctest: +NORMALIZE_WHITESPACE
[array([[0, 6], [0, 0]]), array([[0, 4], [0, 0]]), array([[0, 2], [0, 0]])]
>>> diagrams_generation(2, False, 2, True)  #doctest: +NORMALIZE_WHITESPACE
[array([[0, 2, 2], [0, 0, 2], [0, 0, 0]]),
 array([[0, 1, 1], [0, 0, 3], [0, 0, 0]])]
adg.bmbpt.order_diagrams(diagrams)[source]

Order the BMBPT diagrams and return number of diags for each type.

Parameters:diagrams (list) – Possibly redundant BmbptFeynmanDiagrams.
Returns:
First element is the list of topologically unique, ordered
diagrams. Second element is a dict with the number of diagrams for each major type.
Return type:(tuple)
adg.bmbpt.produce_expressions(diagrams, diagrams_time)[source]

Produce and store the expressions associated to the BMBPT diagrams.

Parameters:
  • diagrams (list) – The list of all BmbptFeynmanDiagrams.
  • diagrams_time (list) – Their associates TSDs.
adg.bmbpt.write_header(tex_file, commands, diags_nbs)[source]

Write overall header for BMBPT result file.

Parameters:
  • tex_file (file) – The ouput LaTeX file of the program.
  • commands (Namespace) – Flags for the program run.
  • diags_nbs (dict) – The number of diagrams per type.