Public API Reference
PowerNetworkMatrices.ABA_Matrix
— TypeStructure containing the ABA matrix and other relevant data.
Arguments
data::SparseArrays.SparseMatrixCSC{Float64, Int}
: the ABA matrix coming from the product between the Incidence Matrix A and the Matrix BA.axes<:NTuple{2, Dict}
: Tuple containing two identical vectors, both containing the number of each bus of the network (each one related to a row/column of the Matrix in "data"), excluding the slack buses.lookup<:NTuple{2, Dict}
: Tuple containing 2 Dictionaries mapping the number of rows and columns with the number of the buses.ref_bus_positions::Set{Int}
: Vector containing the indexes of the columns of the BA matrix corresponding to the reference busesK<:Union{Nothing, KLU.KLUFactorization{Float64, Int}}
: either nothing or a container for KLU factorization matrices (LU factorization)radial_network_reduction::RadialNetworkReduction
: Structure containing the radial branches and leaf buses that were removed while evaluating the matrix
PowerNetworkMatrices.ABA_Matrix
— MethodABA_Matrix(sys; factorize, reduce_radial_branches)
Builds the ABA matrix from a System
Arguments
sys::PSY.System
: system to consider
Keyword arguments
factorize
: if true populates ABA_Matrix.K with KLU factorization matrices
PowerNetworkMatrices.AdjacencyMatrix
— TypeNodal incidence matrix (Adjacency) is an N x N matrix describing a power system with N buses. It represents the directed connectivity of the buses in a power system.
The AdjacencyMatrix Struct is indexed using the Bus Numbers, no need for them to be sequential
Arguments
data::SparseArrays.SparseMatrixCSC{Int8, Int}
: stores the incidence matrixaxes<:NTuple{2, Dict}
: Tuple containing two vectors, the first one contains the names of each line of the network (each one related to a row of the Matrix in "data"), the second one contains the names of each bus of the network (each one related to a column of the Matrix in "data")lookup<:NTuple{2, Dict}
: Tuple containing 2 Dictionaries mapping the number of rows and columns with the names of branches and busesref_bus_positions::Set{Int}
: Vector containing the indexes of the columns of the BA matrix corresponding to the reference buses
PowerNetworkMatrices.AdjacencyMatrix
— MethodAdjacencyMatrix(
branches,
buses;
check_connectivity,
kwargs...
)
Builds a AdjacencyMatrix from a collection of buses and branches. The return is an N x N AdjacencyMatrix Array indexed with the bus numbers.
Arguments
check_connectivity::Bool
: Checks connectivity of the network using Depth First Search (DFS) algorithm
PowerNetworkMatrices.AdjacencyMatrix
— MethodAdjacencyMatrix(sys; check_connectivity, kwargs...)
Builds a AdjacencyMatrix from the system. The return is an N x N AdjacencyMatrix Array indexed with the bus numbers.
Arguments
check_connectivity::Bool
: Checks connectivity of the network using Depth First Search (DFS) algorithm
PowerNetworkMatrices.BA_Matrix
— TypeStructure containing the BA matrix and other relevant data.
Arguments
data::SparseArrays.SparseMatrixCSC{Float64, Int}
: the transposed BA matrix coming from the product between the Incidence Matrix A and the Matrix of Susceptance Baxes<:NTuple{2, Dict}
: Tuple containing two vectors, the first one contains the names of each buse of the network (each one related to a row of the Matrix in "data"), the second one contains the names of each line of the network (each one related to a column of the Matrix in "data")lookup<:NTuple{2, Dict}
: Tuple containing 2 Dictionaries mapping the number of rows and columns with the names of buses and branchesref_bus_positions::Set{Int}
: Set containing the indexes of the columns of the BA matrix corresponding to the reference busesradial_network_reduction::RadialNetworkReduction
: Structure containing the radial branches and leaf buses that were removed while evaluating the matrix
PowerNetworkMatrices.BA_Matrix
— MethodBA_Matrix(sys; reduce_radial_branches)
Build the BA matrix from a given System
Arguments
sys::PSY.System
: PSY system for which the matrix is constructedreduce_radial_branches::Bool
: if True the matrix is build considering radial branches removed from the system
PowerNetworkMatrices.IncidenceMatrix
— TypeIncidence matrix: shows connection between buses, defining lines
Arguments
data::SparseArrays.SparseMatrixCSC{Int8, Int}
: the actual Incidence matrix.axes<:NTuple{2, Dict}
: Tuple containing two vectors (the first one showing the branches names, the second showing the buses numbers).lookup<:NTuple{2, Dict}
: Tuple containing two dictionaries, the first mapping the branches and buses with their enumerated indexes.ref_bus_positions::Set{Int}
: Vector containing the indices of the reference slack buses.radial_network_reduction::RadialNetworkReduction
: Structure containing the radial branches and leaf buses that were removed while evaluating the matrix
PowerNetworkMatrices.IncidenceMatrix
— MethodIncidenceMatrix(sys)
Builds the Incidence matrix of the system by evaluating the actual matrix and other relevant values.
Arguments
sys::PSY.System
: the PowerSystem system to considerreduce_radial_branches::Bool
: if True the matrix will be evaluated discarding all the radial branches and leaf buses (optional, default value is false)
PowerNetworkMatrices.LODF
— TypeThe Line Outage Distribution Factor (LODF) matrix gathers a sensitivity coefficients of how a change in a line's flow affects the flows on other lines in the system.
Arguments
data<:AbstractArray{Float64, 2}
: the transposed LODF matrix.axes<:NTuple{2, Dict}
: Tuple containing two identical vectors containing the names of the branches related to each row/column.lookup<:NTuple{2, Dict}
: Tuple containing two identical dictionaries mapping the branches their enumerated indexes (row and column numbers).tol::Base.RefValue{Float64}
: tolerance used for sparsifying the matrix (dropping element whose absolute value is below this threshold).radial_network_reduction::RadialNetworkReduction
: Structure containing the radial branches and leaf buses that were removed while evaluating the matrix
PowerNetworkMatrices.LODF
— MethodLODF(
branches,
buses;
linear_solver,
tol,
radial_network_reduction
)
Builds the LODF matrix given the data of branches and buses of the system.
Arguments
branches
: vector of the System AC branchesbuses::Vector{PSY.ACBus}
: vector of the System buses
Keyword Arguments
linear_solver
::String linear solver to use for matrix evaluation. Available options: "KLU", "Dense" (OpenBLAS) or "MKLPardiso". Default solver: "KLU".tol::Float64
: Tolerance to eliminate entries in the LODF matrix (default eps())radial_network_reduction::RadialNetworkReduction
: Structure containing the radial branches and leaf buses that were removed while evaluating the ma
PowerNetworkMatrices.LODF
— MethodLODF(A, ABA, BA; linear_solver, tol, reduce_radial_branches)
Builds the LODF matrix given the Incidence Matrix and the PTDF matrix of the system.
NOTE: this method does not support distributed slack bus.
Arguments
A::IncidenceMatrix
: Structure containing the Incidence matrix of the system.ABA::ABA_Matrix
: Structure containing the ABA matrix of the system.BA::BA_Matrix
: Structure containing the transposed BA matrix of the system.linear_solver::String
: Linear solver to be used. Options are "Dense" and "KLU".tol::Float64
: Tolerance to eliminate entries in the LODF matrix (default eps()).reduce_radial_branches::Bool
: True to reduce the network by simplifying the radial branches and mapping the eliminate buses
PowerNetworkMatrices.LODF
— MethodLODF(A, PTDFm; linear_solver, tol, reduce_radial_branches)
Builds the LODF matrix given the Incidence Matrix and the PTDF matrix of the system.
NOTE: tol is referred to the LODF sparsification, not the PTDF one. PTDF matrix must be considered as NON sparsified ("tol" argument not specified when calling the PTDF method).
Arguments
A::IncidenceMatrix
: Structure containing the Incidence matrix of the system.PTDFm::PTDF
: Strucutre containing the transposed PTDF matrix of the system.linear_solver::String
: Linear solver to be used. Options are "Dense" and "KLU".tol::Float64
: Tolerance to eliminate entries in the LODF matrix (default eps()).reduce_radial_branches::Bool
: True to reduce the network by simplifying the radial branches and mapping the eliminate buses
PowerNetworkMatrices.LODF
— MethodLODF(sys; reduce_radial_branches, kwargs...)
Builds the LODF matrix from a system. Note that reduce_radial_branches
kwargs is explicitly mentioned because needed inside of the function.
Arguments
sys::PSY.System
: Power Systems system
Keyword Arguments
reduce_radial_branches::Bool=false
: if True the matrix will be evaluated discarding all the radial branches and leaf buses (optional, default value is false)
PowerNetworkMatrices.PTDF
— TypePower Transfer Distribution Factors (PTDF) indicate the incremental change in real power that occurs on transmission lines due to real power injections changes at the buses.
The PTDF struct is indexed using the Bus numbers and Branch names.
Arguments
data<:AbstractArray{Float64, 2}
: the transposed PTDF matrix.axes<:NTuple{2, Dict}
: Tuple containing two vectors: the first one showing the bus numbers, the second showing the branch names. The information contained in this field matches the axes of the fielsdata
.lookup<:NTuple{2, Dict}
: Tuple containing two dictionaries mapping the bus numbers and branch names with the indices of the matrix contained indata
.subnetworks::Dict{Int, Set{Int}}
: dictionary containing the set of bus indexes defining the subnetworks of the system.tol::Base.RefValue{Float64}
: tolerance used for sparsifying the matrix (dropping element whose absolute value is below this threshold).radial_network_reduction::RadialNetworkReduction
: Structure containing the radial branches and leaf buses that were removed while evaluating the matrix
PowerNetworkMatrices.PTDF
— MethodPTDF(filename)
Deserialize a PTDF from an HDF5 file.
Arguments
filename::AbstractString
: File containing a serialized PTDF.
PowerNetworkMatrices.PTDF
— MethodPTDF(
branches,
buses;
dist_slack,
linear_solver,
tol,
radial_network_reduction
)
Builds the PTDF matrix from a group of branches and buses. The return is a PTDF array indexed with the bus numbers.
Arguments
branches
: vector of the System AC branchesbuses::Vector{PSY.ACBus}
: vector of the System buses
Keyword Arguments
dist_slack::Vector{Float64}
: vector of weights to be used as distributed slack bus. The distributed slack vector has to be the same length as the number of buseslinear_solver::String
: Linear solver to be used. Options are "Dense", "KLU" and "MKLPardisotol::Float64
: Tolerance to eliminate entries in the PTDF matrix (default eps())radial_network_reduction::RadialNetworkReduction
: Structure containing the radial branches and leaf buses that were removed while evaluating the ma
PowerNetworkMatrices.PTDF
— MethodPTDF(
A,
BA;
dist_slack,
linear_solver,
tol,
reduce_radial_branches
)
Builds the PTDF matrix from a system. The return is a PTDF array indexed with the bus numbers.
Arguments
A::IncidenceMatrix
: Incidence Matrix (full structure)BA::BA_Matrix
: BA matrix (full structure)
Keyword Arguments
dist_slack::Vector{Float64}
: Vector of weights to be used as distributed slack bus. The distributed slack vector has to be the same length as the number of buses.linear_solver::String
: Linear solver to be used. Options are "Dense", "KLU" and "MKLPardiso.tol::Float64
: Tolerance to eliminate entries in the PTDF matrix (default eps()).reduce_radial_branches::Bool
: True to reduce the network by simplifying the radial branches and mapping the eliminate buses
PowerNetworkMatrices.PTDF
— MethodPTDF(sys; dist_slack, reduce_radial_branches, kwargs...)
Builds the PTDF matrix from a system. The return is a PTDF array indexed with the bus numbers. Note that dist_slack
and reduce_radial_branches
kwargs are explicitly mentioned because needed inside of the function.
Arguments
sys::PSY.System
: PSY system for which the matrix is constructed
Keyword Arguments
dist_slack::Vector{Float64}=Float64[]
: vector of weights to be used as distributed slack bus. The distributed slack vector has to be the same length as the number of busereduce_radial_branches::Bool=false
: if True the matrix will be evaluated discarding all the radial branches and leaf buses (optional, default value is false)
PowerNetworkMatrices.RadialNetworkReduction
— MethodRadialNetworkReduction(A)
Structure to collect information the branches in the system that are radial and can be ignored in the models.
Arguments
A::IncidenceMatrix
: IncidenceMatrix
PowerNetworkMatrices.VirtualLODF
— TypeThe Virtual Line Outage Distribution Factor (VirtualLODF) structure gathers the rows of the LODF matrix as they are evaluated on-the-go. These rows are evalauted independently, cached in the structure and do not require the computation of the whole matrix (therefore significantly reducing the computational requirements).
The VirtualLODF is initialized with no row stored.
The VirtualLODF struct is indexed using branch names.
Arguments
K::KLU.KLUFactorization{Float64, Int}
: LU factorization matrices of the ABA matrix, evaluated by means of KLU.BA::SparseArrays.SparseMatrixCSC{Float64, Int}
: BA matrix.A::SparseArrays.SparseMatrixCSC{Int8, Int}
: Incidence matrix.inv_PTDF_A_diag::Vector{Float64}
: Vector contiaining the element-wise reciprocal of the diagonal elements coming from multuiplying the PTDF matrix with th Incidence matrixref_bus_positions::Set{Int}
: Vector containing the indexes of the rows of the transposed BA matrix corresponding to the reference buses.axes<:NTuple{2, Dict}
: Tuple containing two vectors showing the branch names.lookup<:NTuple{2, Dict}
: Tuple containing two dictionaries, mapping the branches names the enumerated row indexes indexes.valid_ix::Vector{Int}
: Vector containing the row/columns indices of matrices related the buses which are not slack ones.temp_data::Vector{Float64}
: Temporary vector for internal use.cache::RowCache
: Cache were LODF rows are stored.subnetworks::Dict{Int, Set{Int}}
: Dictionary containing the subsets of buses defining the different subnetwork of the system.tol::Base.RefValue{Float64}
: Tolerance related to scarification and values to drop.radial_network_reduction::RadialNetworkReduction
: Structure containing the radial branches and leaf buses that were removed while evaluating the matrix
PowerNetworkMatrices.VirtualLODF
— MethodVirtualLODF(
branches,
buses;
tol,
max_cache_size,
persistent_lines,
radial_network_reduction
)
Builds the LODF matrix from a group of branches and buses. The return is a VirtualLODF struct with an empty cache.
Arguments
branches
: Vector of the system's AC branches.buses::Vector{PSY.ACBus}
: Vector of the system's buses.
Keyword Arguments
tol::Float64 = eps()
: Tolerance related to sparsification and values to drop.max_cache_size::Int
: max cache size in MiB (inizialized as MAXCACHESIZE_MiB).persistent_lines::Vector{String}
: line to be evaluated as soon as the VirtualPTDF is created (initialized as empty vector of strings).radial_network_reduction::RadialNetworkReduction
: Structure containing the radial branches and leaf buses that were removed while evaluating the matrix
PowerNetworkMatrices.VirtualLODF
— MethodVirtualLODF(sys; reduce_radial_branches, kwargs...)
Builds the Virtual LODF matrix from a system. The return is a VirtualLODF struct with an empty cache.
Arguments
sys::PSY.System
: PSY system for which the matrix is constructed
Keyword Arguments
reduce_radial_branches::Bool=false
: if True the matrix will be evaluated discarding all the radial branches and leaf buses (optional, default value is false)kwargs...
: other keyword arguments used by VirtualPTDF
PowerNetworkMatrices.VirtualPTDF
— TypeThe Virtual Power Transfer Distribution Factor (VirtualPTDF) structure gathers the rows of the PTDF matrix as they are evaluated on-the-go. These rows are evalauted independently, cached in the structure and do not require the computation of the whole matrix (therefore significantly reducing the computational requirements).
The VirtualPTDF is initialized with no row stored.
The VirtualPTDF is indexed using branch names and bus numbers as for the PTDF matrix.
Arguments
K::KLU.KLUFactorization{Float64, Int}
: LU factorization matrices of the ABA matrix, evaluated by means of KLUBA::SparseArrays.SparseMatrixCSC{Float64, Int}
: BA matricref_bus_positions::Set{Int}
: Vector containing the indexes of the columns of the BA matrix corresponding to the reference busesdist_slack::Vector{Float64}
: Vector of weights to be used as distributed slack bus. The distributed slack vector has to be the same length as the number of buses.axes<:NTuple{2, Dict}
: Tuple containing two vectors: the first one showing the branches names, the second showing the buses numbers. There is no link between the order of the vector of the branches names and the way the PTDF rows are stored in the cache.lookup<:NTuple{2, Dict}
: Tuple containing two dictionaries, mapping the branches and buses with their enumerated indexes. The branch indexes refer to the key of the cache dictionary. The bus indexes refer to the position of the elements in the PTDF row stored.temp_data::Vector{Float64}
: Temporary vector for internal use.valid_ix::Vector{Int}
: Vector containing the row/columns indices of matrices related the buses which are not slack ones.cache::RowCache
: Cache were PTDF rows are stored.subnetworks::Dict{Int, Set{Int}}
: Dictionary containing the subsets of buses defining the different subnetwork of the system.tol::Base.RefValue{Float64}
: Tolerance related to scarification and values to drop.radial_network_reduction::RadialNetworkReduction
: Structure containing the radial branches and leaf buses that were removed while evaluating the matrix
PowerNetworkMatrices.VirtualPTDF
— MethodVirtualPTDF(
branches,
buses;
dist_slack,
tol,
max_cache_size,
persistent_lines,
radial_network_reduction
)
Builds the PTDF matrix from a group of branches and buses. The return is a VirtualPTDF struct with an empty cache.
Arguments
branches
: Vector of the system's AC branches.buses::Vector{PSY.ACBus}
: Vector of the system's buses.
Keyword Arguments
dist_slack::Vector{Float64} = Float64[]
: Vector of weights to be used as distributed slack bus. The distributed slack vector has to be the same length as the number of buses.tol::Float64 = eps()
: Tolerance related to sparsification and values to drop.max_cache_size::Int
: max cache size in MiB (inizialized as MAXCACHESIZE_MiB).persistent_lines::Vector{String}
: line to be evaluated as soon as the VirtualPTDF is created (initialized as empty vector of strings).radial_network_reduction::RadialNetworkReduction
: Structure containing the radial branches and leaf buses that were removed while evaluating the matrix
PowerNetworkMatrices.VirtualPTDF
— MethodVirtualPTDF(
sys;
dist_slack,
reduce_radial_branches,
kwargs...
)
Builds the Virtual PTDF matrix from a system. The return is a VirtualPTDF struct with an empty cache.
Arguments
sys::PSY.System
: PSY system for which the matrix is constructed
Keyword Arguments
dist_slack::Vector{Float64}=Float64[]
: vector of weights to be used as distributed slack bus. The distributed slack vector has to be the same length as the number of busereduce_radial_branches::Bool=false
: if True the matrix will be evaluated discarding all the radial branches and leaf buses (optional, default value is false)kwargs...
: other keyword arguments used by VirtualPTDF
PowerNetworkMatrices.Ybus
— TypeYbus(branches, buses; ...)
Ybus(branches, buses, fixed_admittances; check_connectivity)
Builds a Ybus from a collection of buses and branches. The return is a Ybus Array indexed with the bus numbers and the branch names.
Arguments
check_connectivity::Bool
: Checks connectivity of the network using Depth First Search (DFS)
PowerNetworkMatrices.Ybus
— TypeNodal admittance matrix (Ybus) is an N x N matrix describing a power system with N buses. It represents the nodal admittance of the buses in a power system.
The Ybus Struct is indexed using the Bus Numbers, no need for them to be sequential
PowerNetworkMatrices.Ybus
— MethodYbus(sys; kwargs...)
Builds a Ybus from the system. The return is a Ybus Array indexed with the bus numbers and the branch names.
Arguments
check_connectivity::Bool
: Checks connectivity of the network
PowerNetworkMatrices.factorize
— Methodfactorize(ABA)
Evaluates the LU factorization matrices of the ABA matrix, using KLU.
Arguments
ABA::ABA_Matrix{Ax, L, Nothing} where {Ax, L <: NTuple{2, Dict}}
: container for the ABA matrix, with ABA.K == nothing (LU matrices in K not evaluated)
PowerNetworkMatrices.find_subnetworks
— Methodfind_subnetworks(M)
Evaluates subnetworks by looking for the subsets of buses connected each other, but not connected with buses of other subsets.
PowerNetworkMatrices.find_subnetworks
— Methodfind_subnetworks(sys)
Finds the subnetworks in a system using Depth First Search (DFS). Returns a dictionary keyed by the reference bus of the subnetworks if they exist
PowerNetworkMatrices.find_subnetworks
— Methodfind_subnetworks(M, bus_numbers)
Finds the subnetworks present in the considered System. This is evaluated by taking a the ABA or Adjacency Matrix.
Arguments
M::SparseArrays.SparseMatrixCSC
: input sparse matrix.bus_numbers::Vector{Int}
: vector containing the indices of the system's buses.
PowerNetworkMatrices.from_hdf5
— Methodfrom_hdf5(_, filename)
Deserialize a PTDF from an HDF5 file.
Arguments
::Type{PTDF}
:filename::AbstractString
: File containing a serialized PTDF.
PowerNetworkMatrices.to_hdf5
— Methodto_hdf5(ptdf, filename; compress, compression_level, force)
Serialize the PTDF to an HDF5 file.
Arguments
ptdf::PTDF
: matrixfilename::AbstractString
: File to createcompress::Bool
: Whether to enabled compression, defaults to true.compression_level::Int
: Compression level to use if compression is enabled.force::Bool
: Whether to overwrite the file if it already exists, defaults to false.
PowerNetworkMatrices.validate_connectivity
— Methodvalidate_connectivity(M)
Validates connectivity by checking that the number of subnetworks is 1 (fully connected network).
PowerNetworkMatrices.validate_connectivity
— Methodvalidate_connectivity(sys)
Checks the network connectivity of the system using Depth First Search (DFS)