ITensors functions
Utility ITensors functions typically used by the TNDecode module and the TNDistance module.
TensorNetworkCodes.code_to_Itensor — Functioncode_to_Itensor(code::QuantumCode,logical_indices::Array{Index{Int64},1},physical_indices::Array{Index{Int64},1})
-> ITensorReturns an ITensor describing the logical cosets of the code. Indices (Index{Int64}) correspond to logical and physical qubits of the code.
Examples
julia> using ITensors;
julia> logical_indices = [Index(4,"logical")];
julia> physical_indices = [Index(4,"physical") for _ in 1:5];
julia> tensor = code_to_Itensor(five_qubit_surface_code(),logical_indices,physical_indices);
julia> dims(tensor) # tensor has six legs
(4, 4, 4, 4, 4, 4)TensorNetworkCodes.identity_coset — Functionidentity_coset(tensor::ITensor) -> ITensorGiven an ITensor descibing a code, return the ITensor describing only the identity coset, i.e., the stabilizer group.
Examples
julia> using ITensors;
julia> logical_indices = [Index(4,"logical")];
julia> physical_indices = [Index(4,"physical") for _ in 1:7];
julia> tensor = code_to_Itensor(steane_code(),logical_indices,physical_indices);
julia> sum(identity_coset(tensor)) # stabilizer group has 64 elements
64.0TensorNetworkCodes.all_cosets — Functionall_cosets(tensor::ITensor) -> ITensorGiven an ITensor descibing a code, return the ITensor describing all the cosets by summing over all logical indices. Works similarly to identity_coset.
TensorNetworkCodes.physical_tensor — Functionphysical_tensor(index::Index{Int64},error_prob::Float64,pauli::Int64) -> ITensorReturns a one leg ITensor describing the error on one physical qubit. Assumes depolarizing noise as the error model. pauli permutes the values (useful for making decoders to account for pure_errors).
Examples
julia> using ITensors;
julia> index = Index(4,"physical");
julia> pauli = 0;
julia> error_prob = 0.3;
julia> a = array(physical_tensor(index,error_prob,pauli));
julia> println(round.(a,digits=1))
[0.7, 0.1, 0.1, 0.1]
julia> pauli = 2;
julia> a = array(physical_tensor(index,error_prob,pauli));
julia> println(round.(a,digits=1))
[0.1, 0.1, 0.7, 0.1]TensorNetworkCodes.create_virtual_tensor — Functioncreate_virtual_tensor(code::TensorNetworkCode,node::Int64) -> ITensorCreates the ITensor describing the seed code at node. Necessary for, e.g., tensor-network decoding or distance calculation.