Internal API
PowerNetworkMatrices.Adjacency
— 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 Adjacency Struct is indexed using the Bus Numbers, no need for them to be sequential
PowerNetworkMatrices.Adjacency
— MethodAdjacency(branches, nodes; check_connectivity, kwargs...)
Builds a Adjacency from a collection of buses and branches. The return is an N x N Adjacency Array indexed with the bus numbers.
Keyword arguments
check_connectivity::Bool
: Checks connectivity of the network using Goderya's algorithm
PowerNetworkMatrices.Adjacency
— MethodAdjacency(sys; check_connectivity, kwargs...)
Builds a Adjacency from the system. The return is an N x N Adjacency Array indexed with the bus numbers.
Keyword arguments
check_connectivity::Bool
: Checks connectivity of the network using Goderya's algorithmconnectivity_method::Function = goderya_connectivity
: method (goderya_connectivity
ordfs_connectivity
) for connectivity validation
PowerNetworkMatrices.PowerNetworkMatrix
— TypeType PowerNetworkMatrix gathers all the different types of Matrices considered in this package
PowerNetworkMatrices.PowerNetworkMatrixKey
— TypeStructure to store the keys of a power network matrix
Arguments
I<:Tuple
: turple containing the indices of the matrix
PowerNetworkMatrices.PowerNetworkMatrixKeys
— TypeStructure to store the keys of a power network matrix
Arguments
product_iter::Base.Iterators.ProductIterator{T} where T <: Tuple
: iterator of the indices of the network power matrix
PowerNetworkMatrices.RowCache
— TypeStructure used for saving the rows of the Virtual PTDF and LODF matrix.
Arguments
temp_cache::Dict{Int, Union{Vector{Float64}, SparseArrays.SparseVector{Float64}}}
: Dictionay saving the row of the PTDF/LODF matrixpersistent_cache_keys::Set{Int}
: Set listing the rows to keep intemp_cache
max_cache_size::Int
Defines the maximum allowed cache size (rows*row_size)max_num_keys::Int
Defines the maximum number of keys saved (rows of the matrix)
PowerNetworkMatrices.RowCache
— MethodRowCache(max_cache_size, persistent_rows, row_size)
Structure used for saving the rows of the Virtual PTDF and LODF matrix.
Arguments
max_cache_size::Int
Defines the maximum allowed cache size (rows*row_size).persistent_rows::Set{Int}
: Set listing the rows to keep intemp_cache
.row_size
Defines the size of the single row to store.
Base.eachindex
— Methodeachindex(vlodf)
Gives the cartesian indexes of the LODF matrix.
Base.eachindex
— Methodeachindex(vptdf)
Gives the cartesian indexes of the PTDF matrix (same as the BA one).
Base.empty!
— Methodempty!(cache)
Erases the cache.
Base.getindex
— Methodgetindex(cache, key)
Gets the row of the stored matrix in cache.
Arguments
cache::RowCache
: cache where the row vector is going to be savedkey::Int
: row number (corresponding to the enumerated branch index) related to the row vector.
Base.getindex
— Methodgetindex(vlodf, row, column)
Gets the value of the element of the LODF matrix given the row and column indices corresponding to the selected and outage branch respectively. If column
is a Colon then the entire row is returned.
Arguments
vlodf::VirtualLODF
: VirtualLODF struct where to evaluate and store the row values.row
: selected line namecolumn
: outage line name. IfColon
then get the values of the whole row.
Base.getindex
— Methodgetindex(vptdf, row, column)
Gets the value of the element of the PTDF matrix given the row and column indices corresponding to the branch and buses one respectively. If column
is a Colon then the entire row is returned.
Arguments
vptdf::VirtualPTDF
: VirtualPTDF struct where to evaluate and store the row values.row
: Branch index.column
: Bus index. If Colon then get the values of the whole row.
Base.haskey
— Methodhaskey(cache, key)
Checks if key
is present as a key of the dictionary in cache
Arguments
cache::RowCache
: cache where data is stored.key::Int
: row number (corresponds to the enumerated branch index).
Base.isempty
— Methodisempty(cache)
Check if cache is empty.
Base.isempty
— Methodisempty(vlodf)
Checks if the any of the fields of VirtualLODF is empty.
Base.isempty
— Methodisempty(vptdf)
Checks if the any of the fields of VirtualPTDF is empty.
Base.length
— Methodlength(cache)
Shows the number of rows stored in cache
Base.setindex!
— Methodsetindex!(_, _, _)
!!! STILL TO IMPLEMENT !!!
Base.setindex!
— Methodsetindex!(_, _, idx)
!!! STILL TO IMPLEMENT !!!
Base.setindex!
— Methodsetindex!(_, _, _)
!!! STILL TO IMPLEMENT !!!
Base.setindex!
— Methodsetindex!(_, _, idx)
!!! STILL TO IMPLEMENT !!!
Base.setindex!
— Methodsetindex!(cache, val, key)
Allocates vector as row of the matrix saved in cache.
Arguments
cache::RowCache
: cache where the row vector is going to be savedval::Union{Vector{Float64}, SparseArrays.SparseVector{Float64}}
: vector to be savedkey::Int
: row number (corresponding to the enumerated branch index) related to the input row vector
Base.size
— Methodsize(vlodf)
Shows the size of the whole LODF matrix, not the number of rows stored.
Base.size
— Methodsize(vptdf)
Gives the size of the whole PTDF matrix, not the number of rows stored.
PowerNetworkMatrices._calculate_PTDF_matrix_DENSE
— Method_calculate_PTDF_matrix_DENSE(
A,
BA,
ref_bus_positions,
dist_slack
)
Funciton for internal use only.
Computes the PTDF matrix by means of the LAPACK and BLAS functions for dense matrices.
Arguments
A::Matrix{Int8}
: Incidence MatrixBA::Matrix{T} where {T <: Union{Float32, Float64}}
: BA matrixref_bus_positions::Set{Int}
: vector containing the indexes of the reference slack buses.dist_slack::Vector{Float64})
: vector containing the weights for the distributed slacks.
PowerNetworkMatrices._calculate_PTDF_matrix_KLU
— Method_calculate_PTDF_matrix_KLU(
A,
BA,
ref_bus_positions,
dist_slack
)
Funciton for internal use only.
Computes the PTDF matrix by means of the KLU.LU factorization for sparse matrices.
Arguments
A::SparseArrays.SparseMatrixCSC{Int8, Int}
: Incidence MatrixBA::SparseArrays.SparseMatrixCSC{Float64, Int}
: BA matrixref_bus_positions::Set{Int}
: vector containing the indexes of the reference slack buses.dist_slack::Vector{Float64}
: vector containing the weights for the distributed slacks.
PowerNetworkMatrices._calculate_PTDF_matrix_MKLPardiso
— Method_calculate_PTDF_matrix_MKLPardiso(
A,
BA,
ref_bus_positions,
dist_slack
)
Funciton for internal use only.
Computes the PTDF matrix by means of the MKL Pardiso for dense matrices.
Arguments
A::SparseArrays.SparseMatrixCSC{Int8, Int}
: Incidence MatrixBA::SparseArrays.SparseMatrixCSC{Float64, Int}
: BA matrixref_bus_positions::Set{Int}
: vector containing the indexes of the referece slack buses.dist_slack::Vector{Float64}
: vector containing the weights for the distributed slacks.
PowerNetworkMatrices.assign_reference_buses!
— Methodassign_reference_buses!(subnetworks, ref_buses)
Takes the reference bus numbers and re-assigns the keys in the subnetwork dictionaries to use the reference bus withing each subnetwork.
PowerNetworkMatrices.calculate_ABA_matrix
— Methodcalculate_ABA_matrix(A, BA, ref_bus_positions)
Evaluates the ABA matrix given the System's Incidence matrix (A), BA matrix and reference bus positions.
Arguments
A::SparseArrays.SparseMatrixCSC{Int8, Int}
: Incidence matrix.BA::SparseArrays.SparseMatrixCSC{Float64, Int}
BA matrix.
NOTE:
- evaluates A with "calculateAmatrix", or extract A.data (if A::IncidenceMatrix)
- evaluates BA with "calculateBAmatrix", or extract BA.data (if A::BA_Matrix)
PowerNetworkMatrices.calculate_A_matrix
— Methodcalculate_A_matrix(branches, buses)
Evaluates the Incidence matrix A given the branches and node of a System.
Arguments
branches
: vector containing the branches of the considered system (should be AC branches).buses::Vector{PSY.ACBus}
: vector containing the buses of the considered system.
NOTE:
- the matrix features all the columns, including the ones related to the reference buses (each column is related to a system's bus).
PowerNetworkMatrices.calculate_BA_matrix
— Methodcalculate_BA_matrix(branches, bus_lookup)
Evaluates the transposed BA matrix given the System's banches, reference bus positions and bus_lookup.
Arguments
branches
: vector containing the branches of the considered system (should be AC branches).ref_bus_positions::Set{Int}
: Vector containing the indexes of the columns of the BA matrix corresponding to the reference busesbus_lookup::Dict{Int, Int}
: dictionary mapping the bus numbers with their enumerated indexes.
PowerNetworkMatrices.calculate_PTDF_matrix_DENSE
— Methodcalculate_PTDF_matrix_DENSE(
branches,
buses,
bus_lookup,
dist_slack
)
Computes the PTDF matrix by means of the LAPACK and BLAS functions for dense matrices.
Arguments
branches
: vector of the System AC branchesbuses::Vector{PSY.ACBus}
: vector of the System busesbus_lookup::Dict{Int, Int}
: dictionary mapping the bus numbers with their enumerated indexes.dist_slack::Vector{Float64}
: vector containing the weights for the distributed slacks.
PowerNetworkMatrices.calculate_PTDF_matrix_KLU
— Methodcalculate_PTDF_matrix_KLU(
branches,
buses,
bus_lookup,
dist_slack
)
Computes the PTDF matrix by means of the KLU.LU factorization for sparse matrices.
Arguments
branches
: vector of the System AC branchesbuses::Vector{PSY.ACBus}
: vector of the System busesbus_lookup::Dict{Int, Int}
: dictionary mapping the bus numbers with their enumerated indexes.dist_slack::Vector{Float64}
: vector containing the weights for the distributed slacks.
PowerNetworkMatrices.calculate_PTDF_matrix_MKLPardiso
— Methodcalculate_PTDF_matrix_MKLPardiso(
branches,
buses,
bus_lookup,
dist_slack
)
Computes the PTDF matrix by means of the MKL Pardiso for dense matrices.
Arguments
branches
: vector of the System AC branchesbuses::Vector{PSY.ACBus}
: vector of the System busesbus_lookup::Dict{Int, Int}
: dictionary mapping the bus numbers with their enumerated indexes.dist_slack::Vector{Float64}
: vector containing the weights for the distributed slacks.
PowerNetworkMatrices.calculate_adjacency
— Methodcalculate_adjacency(branches, buses, bus_lookup)
Evaluates the Adjacency matrix given the System's banches, buses and bus_lookup.
NOTE:
- bus_lookup is a dictionary mapping the bus numbers (as shown in the Systems) with their enumerated indxes.
PowerNetworkMatrices.calculate_adjacency
— Methodcalculate_adjacency(branches, buses)
Evaluates the Adjacency matrix given the banches and buses of a given System.
Arguments
branches
: vector containing the branches of the considered system (should be AC branches).buses::Vector{PSY.ACBus}
: vector containing the buses of the considered system.
PowerNetworkMatrices.calculate_radial_branches
— Methodcalculate_radial_branches(
A,
line_map,
bus_map,
ref_bus_positions
)
used to calculate the branches in the system that are radial and can be ignored in the models by exploring the structure of the incidence matrix
Arguments
A::SparseArrays.SparseMatrixCSC{Int8, Int64}
: Data from the IncidenceMatrixline_map::Dict{String, Int}
: Map of Line Name to Matrix Indexbus_map::Dict{Int, Int}
: Map of Bus Name to Matrix Indexref_bus_positions::Set{Int}
: Set containing the indexes of the columns of the BA matrix corresponding to the reference buses
PowerNetworkMatrices.check_cache_size!
— Methodcheck_cache_size!(cache; new_add)
Check saved rows in cache and delete one not belonging to persistent_cache_keys
.
PowerNetworkMatrices.evaluate_A_matrix_values
— Methodevaluate_A_matrix_values(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 considerradial_network_reduction::RadialNetworkReduction
: Structure containing the radial branches and leaf buses that were removed while evaluating the matrix
PowerNetworkMatrices.find_slack_positions
— Methodfind_slack_positions(buses)
Gets the indices of the reference (slack) buses. NOTE:
- the indices corresponds to the columns of zeros belonging to the PTDF matrix.
- BA and ABA matrix miss the columns related to the reference buses.
PowerNetworkMatrices.get_ac_branches
— Functionget_ac_branches(sys)
get_ac_branches(sys, radial_branches)
Gets the AC branches from a given Systems.
PowerNetworkMatrices.get_bus_indices
— Methodget_bus_indices(branch, bus_lookup)
Evaluates the bus indices for the given branch.
Arguments
branch
: system's branchbus_lookup
: dictionary mapping the system's buses and branches
PowerNetworkMatrices.get_buses
— Functionget_buses(sys)
get_buses(sys, bus_reduction_map)
Gets the non-isolated buses from a given System
PowerNetworkMatrices.get_data
— Methodget_data(mat)
returns the raw array data of the PowerNetworkMatrix
PowerNetworkMatrices.get_lookup
— Methodget_lookup(mat)
returns the lookup tuple of the `PowerNetworkMatrix`. The first entry corresponds
to the first dimension and the second entry corresponds to the second dimension. For
instance in Ybus the first dimension is buses and second dimension is buses too, and in
PTDF the first dimension is branches and the second dimension is buses
PowerNetworkMatrices.get_mapped_bus_number
— Methodget_mapped_bus_number(rb, bus_number)
Interface to obtain the parent bus number of a reduced bus when radial branches are eliminated
Arguments
rb::RadialNetworkReduction
: RadialNetworkReduction objectbus_number::Int
: Bus number of the reduced bus
PowerNetworkMatrices.get_mapped_bus_number
— Methodget_mapped_bus_number(rb, bus)
Interface to obtain the parent bus number of a reduced bus when radial branches are eliminated
Arguments
rb::RadialNetworkReduction
: RadialNetworkReduction objectbus::ACBus
: Reduced bus
PowerNetworkMatrices.get_tol
— Methodget_tol(mat)
Gets the tolerance used for sparsifying the rows of the VirtualLODF matrix
PowerNetworkMatrices.get_tol
— Methodget_tol(vptdf)
Gets the tolerance used for sparsifying the rows of the VirtualPTDF matrix
PowerNetworkMatrices.lookup_index
— Methodlookup_index(i, lookup)
Gets bus indices to a certain branch index
PowerNetworkMatrices.lookup_index
— Methodlookup_index(i, lookup)
Gets bus indices to a certain branch name
Arguments
i::PSY.ACBranch
: Power System AC branchlookup::Dict
: Dictionary mapping branch and buses
PowerNetworkMatrices.lookup_index
— Methodlookup_index(i, lookup)
Gets bus indices to a certain branch name
Arguments
i::PSY.ACBranch
: Power System AC branchlookup::Dict
: Dictionary mapping branches and buses
PowerNetworkMatrices.make_ax_ref
— Methodmake_ax_ref(buses)
Evaluates the map linking the system's buses and branches.
Arguments
buses::AbstractVector{PSY.ACBus}
: system's buses
PowerNetworkMatrices.make_ax_ref
— Methodmake_ax_ref(ax)
Checkes if repetitions are present in the dictionary mapping buses and branches.
Arguments
ax::AbstractVector
: generic abstract vector
PowerNetworkMatrices.purge_one!
— Methodpurge_one!(cache)
Deletes a row from the stored matrix in cache not belonging to the persistentcachekeys set.
PowerNetworkMatrices.sparsify
— Methodsparsify(dense_array, tol)
Return a sparse matrix given a dense one by dropping element whose absolute value is above a certain tolerance.
Arguments
- dense_array::Matrix{Float64}`: input matrix (e.g., PTDF matrix).
tol::Float64
: tolerance.
PowerNetworkMatrices.sparsify
— Methodsparsify(dense_array, tol)
Return a sparse vector given a dense one by dropping element whose absolute value is above a certain tolerance.
Arguments
- dense_array::Vector{Float64}`: input vector (e.g., PTDF row from VirtualPTDF).
tol::Float64
: tolerance.
PowerNetworkMatrices.to_index
— Methodto_index(A, idx)
Given the indices, gets the values of the power network matrix considered
PowerNetworkMatrices.validate_linear_solver
— Methodvalidate_linear_solver(linear_solver)
Validates if the selected linear solver is supported.