Internal Function Documentation

Ops

ITensorMPOConstruction.is_fermionicFunction
is_fermionic(ops::NTuple{N,OpID}, n::Int, op_cache_vec::OpCacheVec) where {N} -> Bool

Compute the fermionic parity of the part of ops acting on sites >= n.

The result is true when an odd number of the operators in ops are marked fermionic in op_cache_vec, and false otherwise.

source
ITensorMPOConstruction.are_equalFunction
are_equal(op1::NTuple{N,OpID}, op2::NTuple{N,OpID}, n::Int) where {N} -> Bool

Check whether two products of OpID are equal from site n onward. The tuples must be sorted by decreasing site.

The comparison proceeds from left to right. If both tuples have already moved to sites < n at the same position, the remaining entries are ignored and the tuples are considered equal. Otherwise, entries must match exactly.

source
ITensorMPOConstruction.terms_eq_fromFunction
terms_eq_from(n::Int) -> Function

Return a comparator function that checks whether two products of OpID are equal from site n onward.

The returned function expects tuples sorted by decreasing site and applies the same comparison rule as are_equal: once both tuples have moved to sites < n at the same position, the remaining entries are ignored.

source
ITensorMPOConstruction.sort_fermion_perm!Function
sort_fermion_perm!(ops::AbstractVector{<:OpID}, op_cache_vec::OpCacheVec) -> Int

Sort ops in place by site index and return the fermionic sign generated by the permutation. The sort is stable.

Whenever two fermionic operators are swapped past one another, the returned sign is multiplied by -1. Non-fermionic swaps do not affect the sign.

source
ITensorMPOConstruction.determine_val_typeFunction
determine_val_type(os::OpIDSum)

Return the numeric element type needed to represent the coefficients and cached operator matrices in os.

This is ComplexF64 if either the stored coefficients or any cached operator matrix has nonzero imaginary part; otherwise it is Float64.

source
ITensorMPOConstruction.for_equal_sitesFunction
for_equal_sites(f::Function, ops::AbstractVector{<:OpID}) -> Nothing

Call f(i, j) for each maximal contiguous block of operators in ops that act on the same site, where i and j denote the start and end of the contiguous block.

This assumes ops are already sorted by site, as happens after sorting terms inside OpIDSum.

source
ITensorMPOConstruction.rewrite_in_operator_basis!Function
rewrite_in_operator_basis!(os::OpIDSum, basis_op_cache_vec::OpCacheVec)

Rewrite each same-site operator product in os into a single operator drawn from basis_op_cache_vec.

For each contiguous block of operators acting on the same site, the product matrix is matched against the supplied basis up to an overall scalar factor. The term coefficient is updated by that factor, the block is replaced by the matched basis operator, and the remaining slots are zeroed. An error is thrown if a product cannot be represented by a single basis operator.

source
ITensorMPOConstruction.op_sum_to_opID_sumFunction
op_sum_to_opID_sum(os::OpSum{C}, sites::Vector{<:Index}) -> OpIDSum

Convert an ITensor OpSum into an OpIDSum.

Each distinct local operator encountered on a site is assigned an integer id and cached in an OpCacheVec, with id 1 reserved for the identity. Terms are then added to the resulting OpIDSum using those cached identifiers.

source
ITensorMPOConstruction.check_os_for_errorsFunction
check_os_for_errors(os::OpIDSum) -> Nothing

Validate structural assumptions for os.

This checks that, within every term, operators are sorted by site, at most one operator acts on each site, all terms carry the same total QN flux, and every term has even fermion parity.

source
ITensorMPOConstruction.prepare_opID_sum!Function
prepare_opID_sum!(os::OpIDSum, basis_op_cache_vec::Union{Nothing,OpCacheVec})

Apply optional preprocessing to os before MPO construction.

Currently this rewrites terms into basis_op_cache_vec when a basis cache is provided and otherwise leaves os unchanged.

source
ITensorMPOConstruction.get_onsite_opFunction
get_onsite_op(ops::NTuple{N,OpID{Ti}}, n::Int) where {N, Ti} -> Ti

Return the operator id in ops acting on site n.

If no operator acts on site n, this returns 1 a.k.a the identity operator.

source

Bipartite Graph

ITensorMPOConstruction.BipartiteGraphType
BipartiteGraph{L,R,C}

Weighted bipartite graph with typed left and right vertices.

Type Parameters

  • L: The type of the left vertices.
  • R: The type of the right vertices.
  • C: The scalar edge weight type.

Fields:

  • left_vertices: metadata stored for each left vertex.
  • right_vertices: metadata stored for each right vertex.
  • right_vertex_ids_from_left: adjacency list from each left vertex to the connected right-vertex ids.
  • edge_weights_from_left: edge weights stored in parallel with right_vertex_ids_from_left.
source
ITensorMPOConstruction.BipartiteGraphConnectedComponentsType
BipartiteGraphConnectedComponents

Struct containing the connected components of a BipartiteGraph.

Fields:

  • lvs_of_component: for each connected component, the global ids of its left vertices.
  • component_position_of_rvs: for each global right vertex id, its local id (position) within its connected component, or an unused sentinel for isolated right vertices.
  • rv_size_of_component: number of right vertices in each component.
source
ITensorMPOConstruction.left_sizeFunction
left_size(g::BipartiteGraph) -> Int

Return the number of left vertices in g.

source
left_size(ccs::BipartiteGraphConnectedComponents, cc::Int)

Return the number of left vertices in connected component cc.

source
ITensorMPOConstruction.clear_edges_from_left!Function
clear_edges_from_left!(g::BipartiteGraph, lv_id::Integer) -> Nothing

Remove all edges incident on the left vertex lv_id.

This clears both the adjacent right-vertex ids and their corresponding edge weights, and releases any retained capacity in those per-vertex adjacency lists.

source
ITensorMPOConstruction.combine_duplicate_adjacent_right_vertices!Function
combine_duplicate_adjacent_right_vertices!(g::BipartiteGraph, eq::Function) -> Vector{Int}

Combine adjacent right vertices in g that compare equal under eq.

The right vertices are assumed to already be grouped so that duplicate vertices appear contiguously. The first vertex of each equal run is kept, later vertices in the run are removed, and every right-vertex id stored in the left adjacency lists is remapped to the surviving vertex id. The returned vector maps each original right-vertex id to its new position after compaction.

source
ITensorMPOConstruction.compute_connected_componentsFunction
compute_connected_components(g::BipartiteGraph) -> BipartiteGraphConnectedComponents

Compute the connected components of g, keeping only components which connect the left and right sides of the bipartite graph (i.e. discarding isolated vertices).

The returned object records, for each retained component, which left vertices it contains and how global right-vertex ids map into local positions within their component.

source
ITensorMPOConstruction.get_cc_matrixFunction
get_cc_matrix(g::BipartiteGraph, ccs::BipartiteGraphConnectedComponents, cc::Int; clear_edges=false)
    -> Tuple{SparseMatrixCSC,Vector{Int},Vector{Int}}

Extract the subgraph corresponding to connected component cc as a sparse matrix.

The returned tuple contains:

  • the sparse weighted adjacency matrix of the component, with local left/right numbering,
  • the map from local row indices back to global left-vertex ids,
  • the map from local column indices back to global right-vertex ids.

If clear_edges=true, the consumed adjacency lists are emptied from g as the matrix is assembled.

source

Others

ITensorMPOConstruction.BlockSparseMatrixType
BlockSparseMatrix{C}

Dictionary-backed block-sparse matrix representation used for intermediate MPO tensor storage.

Keys are (left_link, right_link) pairs and values are dense local operator matrices for that block.

source
ITensorMPOConstruction.to_sparse_itensorFunction
to_sparse_itensor(offsets, block_sparse_matrices, inds...; tol=0.0, checkflux=true) -> ITensor

Assemble an ITensor from component-local block-sparse MPO data.

offsets gives the starting outgoing-link offset for each connected component, and block_sparse_matrices stores the dense local blocks keyed by link-index pairs. Entries with magnitude at most tol are dropped during the copy.

source
ITensorMPOConstruction.add_to_local_matrix!Function
add_to_local_matrix!(a, weight, local_op, needs_JW_string) -> Nothing

Accumulate a weighted local operator matrix into a.

If needs_JW_string is true, the contribution is multiplied by the diagonal Jordan-Wigner sign pattern expected for 2-state or 4-state fermionic sites.

source
ITensorMPOConstruction.CoSorterElementType
CoSorterElement{T1,T2}

Pair-like element used by CoSorter so that sorting one array carries along a second array with the same permutation.

Taken from https://discourse.julialang.org/t/how-to-sort-two-or-more-lists-at-once/12073/13

source
ITensorMPOConstruction.CoSorterType
CoSorter{T1,T2,S,C}

Lightweight view that exposes two arrays as a single sortable collection, ordering by sortarray and applying the same swaps to coarray.

source
ITensorMPOConstruction.sparse_qrFunction
sparse_qr(A::SparseMatrixCSC, tol::Real, absolute_tol::Bool)
    -> Tuple{Q,R,prow,pcol,rank}

Compute a sparse QR factorization of A using SuiteSparse SPQR.

If absolute_tol is false, tol is interpreted relative to SPQR's default tolerance scale. The return values are the orthogonal factor Q, upper-triangular factor R, row and column permutations, and the numerical rank.

source
ITensorMPOConstruction.for_non_zeros_batchFunction
for_non_zeros_batch(f::Function, A::SparseMatrixCSC, max_col::Int) -> Nothing

Iterate over the nonzero entries of the first max_col columns of A, calling f(values, rows, col) once per nonempty column.

values and rows are views into the storage of A for that column.

source
for_non_zeros_batch(f::Function, Q::SparseArrays.SPQR.QRSparseQ, max_col::Int) -> Nothing

Iterate over the first max_col columns of a sparse SPQR Q factor, forming each column explicitly and calling f(column, col).

source
ITensorMPOConstruction.at_site!Function
at_site!(ValType, g, n, sites, tol, absolute_tol, op_cache_vec; combine_qn_sectors, output_level=0)
    -> Tuple{MPOGraph,Vector{Int},Vector{BlockSparseMatrix{ValType}},Index}

Process one site of the MPO construction algorithm.

For each connected component of g, this routine forms the sparse adjacency matrix, computes a sparse QR factorization, assembles the local MPO tensor blocks for site n, and builds the graph passed to the next site. The connected components are iterated over using threads. The returned tuple contains:

  • the graph for site n + 1,
  • offsets locating each connected component inside the outgoing bond space,
  • the block-sparse local MPO tensors for each component,
  • the outgoing link Index, optionally grouped into merged QN sectors.
source
ITensorMPOConstruction.MPOGraphType
MPOGraph{N,C,Ti}

Type alias for the bipartite graph representation used during MPO construction.

Left vertices store LeftVertex metadata, right vertices store fixed-width tuples of OpIDs describing the remaining operator content of a term, and edge weights carry the scalar coefficients.

source
ITensorMPOConstruction.LeftVertexType
LeftVertex

Represents a left vertex of the bipartite MPO graph. Fields:

  • link: integer label for the incoming (left) MPO link.
  • op_id: identifier specifying the operator acting on the current site.
  • needs_JW_string: whether a Jordan-Wigner string must be inserted when connecting through this vertex.
source
ITensorMPOConstruction.add_to_next_graph!Function
add_to_next_graph!(next_graph, cur_graph, op_cache_vec, cur_site, cur_offset, next_edges) -> Nothing

Append the left vertices and adjacency lists described by next_edges to next_graph.

Each nonempty (bond_index, op_id) entry in next_edges creates one LeftVertex, reusing that entry's (right_vertex_ids, edge_weights) vectors as the outgoing adjacency list. The stored needs_JW_string flag is inferred from the fermionic parity of the connected right vertices.

source
ITensorMPOConstruction.find_first_eq_rvFunction
find_first_eq_rv(g::MPOGraph, j::Int, n::Int) -> Int

Walk backward from right-vertex index j to the first earlier right vertex whose operator tuple is equivalent from site n onward.

This is used to "merge" equivalent right vertices after peeling off the operator acting on the current site.

source
ITensorMPOConstruction.merge_qn_sectorsFunction
merge_qn_sectors(qi_of_cc::Vector{Pair{QN,Int}})
    -> Tuple{Vector{Int},Vector{Pair{QN,Int}}}

Sort connected components by QN sector and merge adjacent entries with the same QN by summing their dimensions.

The first returned vector is the permutation that reorders components, and the second is the merged QN-sector description.

source
ITensorMPOConstruction.build_next_edges_specialization!Function
build_next_edges_specialization!(next_edges, g, cur_site, right_vertex_ids, edge_weights) -> Nothing

Fast path for building outgoing edges when a connected component has only a single left vertex.

For each current edge, this extracts the operator acting on cur_site + 1, finds the unique right vertex from cur_site + 2 onward, and stores the resulting id/weight entry in next_edges.

source
ITensorMPOConstruction.process_single_left_vertex_cc!Function
process_single_left_vertex_cc!(
  matrix_of_cc, qi_of_cc, rank_of_cc, next_edges_of_cc, g, ccs, cc, n, sites,
  has_qns, op_cache_vec
) -> Nothing

Handle the common connected-component case with exactly one left vertex.

This fills the local MPO tensor contribution for the component, records its QN sector and rank, and either applies the terminal scaling on the last site or builds the outgoing edges for the next site.

source
ITensorMPOConstruction.@time_ifMacro
@time_if output_level min_output_level msg ex

Conditionally time and print the execution of ex depending on the current verbosity level.

If output_level > min_output_level, this expands to a timed evaluation using @time, prefixing the printed message with min_output_level levels of two-space indentation followed by msg. Otherwise, ex is evaluated without timing output.

source