Public API Reference
PowerNetworkMatrices.ABA_Matrix — Type
Structure containing the ABA matrix and related power system analysis data.
The ABA matrix represents the bus susceptance matrix computed as A^T * B * A, where A is the incidence matrix and B is the branch susceptance matrix. This matrix is fundamental for DC power flow analysis, sensitivity calculations, and linear power system studies.
Fields
data::SparseArrays.SparseMatrixCSC{Float64, Int}: The ABA matrix data representing the bus susceptance matrix. This square matrix has dimensions equal to the number of buses excluding reference busesaxes::Ax: Tuple containing identical bus number vectors for rows and columns, excluding reference buseslookup::L <: NTuple{2, Dict}: Tuple of identical dictionaries providing fast lookup from bus numbers to matrix indicessubnetwork_axes::Dict{Int, Ax}: Mapping from reference bus numbers to their corresponding subnetwork axesref_bus_position::Vector{Int}: Vector containing the original indices of reference buses before matrix reductionK::F <: Union{Nothing, KLU.KLUFactorization{Float64, Int}}: Optional KLU factorization object for efficient linear system solving. Nothing if unfactorizednetwork_reduction_data::NetworkReductionData: Container for network reduction information applied during matrix construction
Mathematical Properties
- Matrix Form: ABA = A^T * B * A (bus susceptance matrix)
- Dimensions: (nbuses - nref) × (nbuses - nref)
- Symmetry: Positive definite symmetric matrix (for connected networks)
- Sparsity: Inherits sparsity pattern from network topology
Notes
- Reference buses are excluded from the matrix to ensure invertibility
- Factorization enables efficient solving of linear systems Ax = b
- Used primarily for DC power flow analysis and power system sensitivity studies
- Supports various network reduction techniques for computational efficiency
PowerNetworkMatrices.ABA_Matrix — Method
ABA_Matrix(sys; factorize, network_reductions, kwargs...)
ABA_Matrix(sys::PSY.System; factorize::Bool = false, network_reductions::Vector{NetworkReduction} = NetworkReduction[], kwargs...)Construct an ABA_Matrix from a PowerSystems.System by computing A^T * B * A where A is the incidence matrix and B is the branch susceptance matrix. The resulting matrix is fundamental for DC power flow analysis and power system sensitivity studies.
Arguments
sys::PSY.System: The power system from which to construct the ABA matrix
Keyword Arguments
factorize::Bool = false: Whether to perform KLU factorization during construction for efficient linear system solvingnetwork_reductions::Vector{NetworkReduction} = NetworkReduction[]: Vector of network reduction algorithms to apply before matrix constructioninclude_constant_impedance_loads::Bool=true: Whether to include constant impedance loads as shunt admittances in the network modelsubnetwork_algorithm=iterative_union_find: Algorithm used for identifying electrical islands and connected components- Additional keyword arguments are passed to the underlying
Ybusconstructor
Returns
ABA_Matrix: The constructed ABA matrix structure containing:- Bus susceptance matrix data (excluding reference buses)
- Network topology information and reference bus positions
- Optional KLU factorization for efficient solving
Mathematical Process
- Ybus Construction: Creates admittance matrix from system data
- Incidence Matrix: Computes bus-branch incidence matrix A
- BA Matrix: Forms branch susceptance weighted incidence matrix
- ABA Computation: Calculates A^T * B * A (bus susceptance matrix)
- Reference Bus Removal: Excludes reference buses for invertibility
- Optional Factorization: Performs KLU decomposition if requested
Notes
- Reference buses are automatically detected and excluded from the final matrix
- Factorization significantly improves performance for repeated linear system solves
- Network reductions can dramatically improve computational efficiency for large systems
- The resulting matrix supports PTDF, LODF, and other power system analysis calculations
PowerNetworkMatrices.AdjacencyMatrix — Type
AdjacencyMatrix{Ax, L} <: PowerNetworkMatrix{Int8}An N × N adjacency matrix representing the connectivity structure of a power system with N buses. This matrix describes the directed connectivity between buses, where non-zero entries indicate electrical connections through transmission lines, transformers, or other network elements.
The matrix is indexed using bus numbers, which do not need to be sequential. Each element A[i,j] is non-zero if there is a direct electrical connection between bus i and bus j. Diagonal elements are typically zero since self-loops are not meaningful in power network topology.
Fields
data::SparseArrays.SparseMatrixCSC{Int8, Int}: The sparse adjacency matrix storing connectivity information as Int8 values (zero for no connection, non-zero for connection)axes::Ax: Tuple containing the axis labels for both dimensions. The first element contains bus identifiers for rows, the second contains bus identifiers for columns (typically identical)lookup::L: Tuple of dictionaries providing bidirectional mapping between bus numbers and their corresponding matrix indicessubnetwork_axes::Dict{Int, Ax}: Dictionary mapping subnetwork identifiers to their corresponding axis information, used for handling electrical islandsnetwork_reduction_data::NetworkReductionData: Container for network reduction algorithms and their associated data, enabling efficient matrix operations on reduced networks
Examples
# Create from a PowerSystems.System
sys = System("case5.m")
adj = AdjacencyMatrix(sys)
# Create from a Ybus matrix
ybus = Ybus(sys)
adj = AdjacencyMatrix(ybus)
# Check connectivity
is_connected = validate_connectivity(adj)
subnetworks = find_subnetworks(adj)See also: Ybus, IncidenceMatrix, PowerNetworkMatrix
PowerNetworkMatrices.AdjacencyMatrix — Method
AdjacencyMatrix(sys; kwargs...)
AdjacencyMatrix(sys::PSY.System; kwargs...)Construct an AdjacencyMatrix from a PowerSystems.System.
Arguments
sys::PSY.System: The power system from which to construct the adjacency matrix
Keyword arguments
network_reductions::Vector{NetworkReduction}=[]: Network reduction algorithms to applyinclude_constant_impedance_loads::Bool=true: Whether to include constant impedance loads as shunt admittancessubnetwork_algorithm=iterative_union_find: Algorithm for finding electrical islands
Returns
AdjacencyMatrix: An N x N adjacency matrix indexed with bus numbers showing connectivity
PowerNetworkMatrices.AdjacencyMatrix — Method
AdjacencyMatrix(ybus)
AdjacencyMatrix(ybus::Ybus)Construct an AdjacencyMatrix from a Ybus matrix.
Arguments
ybus::Ybus: The Ybus matrix from which to construct the adjacency matrix
Returns
AdjacencyMatrix: The constructed adjacency matrix showing bus connectivity
PowerNetworkMatrices.ArcAdmittanceMatrix — Type
Arc admittance matrix
Arguments
data::SparseArrays.SparseMatrixCSC{ComplexF32, Int}: The arc admittance matrix in the from-to directionaxes<:NTuple{2, Dict}: Tuple containing two vectors (the first one showing the arc tuples, the second showing the buses numbers).lookup<:NTuple{2, Dict}: Tuple containing two dictionaries, the first mapping the arc tuples and the second the buses with their enumerated indexes.network_reduction::NetworkReduction: Structure containing the details of the network reduction applied when computing the matrixdirection::Symbol: Direction of admittance (:FromTo or :ToFrom)
PowerNetworkMatrices.BA_Matrix — Type
Structure containing the BA matrix and related network topology data.
The BA matrix represents the branch-bus incidence matrix weighted by branch susceptances, computed as the product of the incidence matrix A and the susceptance matrix B.
Fields
data::SparseArrays.SparseMatrixCSC{Float64, Int}: The transposed BA matrix data. Each row corresponds to a bus and each column corresponds to a branch, with values representing weighted branch susceptancesaxes::Ax: Tuple containing two vectors: bus numbers (rows) and branch identifiers (columns)lookup::L <: NTuple{2, Dict}: Tuple of dictionaries providing fast lookup from bus/branch names to matrix indicessubnetwork_axes::Dict{Int, Ax}: Mapping from reference bus numbers to their corresponding subnetwork axesnetwork_reduction_data::NetworkReductionData: Container for network reduction information applied during matrix construction
Notes
- The matrix is stored in transposed form for computational efficiency
- Reference buses are identified through
subnetwork_axeskeys - Supports various network reduction techniques for computational efficiency
PowerNetworkMatrices.BA_Matrix — Method
BA_Matrix(sys; network_reductions, kwargs...)
BA_Matrix(sys::PSY.System; network_reductions::Vector{NetworkReduction} = Vector{NetworkReduction}(), kwargs...)Construct a BA_Matrix from a PowerSystems.System by first building the underlying Ybus matrix and then computing the branch-bus incidence matrix weighted by branch susceptances.
Arguments
sys::PSY.System: The power system from which to construct the BA matrix
Keyword Arguments
network_reductions::Vector{NetworkReduction} = Vector{NetworkReduction}(): Vector of network reduction algorithms to apply before matrix constructioninclude_constant_impedance_loads::Bool=true: Whether to include constant impedance loads as shunt admittances in the network modelsubnetwork_algorithm=iterative_union_find: Algorithm used for identifying electrical islands and connected components- Additional keyword arguments are passed to the underlying
Ybusconstructor
Returns
BA_Matrix: The constructed BA matrix structure containing the transposed branch-bus incidence matrix weighted by susceptances, along with network topology information
Notes
- This constructor creates a
Ybusmatrix internally and then converts it to aBA_Matrix - Network reductions can significantly improve computational efficiency for large systems
- The resulting matrix supports DC power flow calculations and sensitivity analysis
PowerNetworkMatrices.BA_Matrix — Method
BA_Matrix(ybus)
BA_Matrix(ybus::Ybus)Construct a BA_Matrix from a Ybus matrix.
Arguments
ybus::Ybus: The Ybus matrix from which to construct the BA matrix
Returns
BA_Matrix: The constructed BA matrix structure containing the transposed BA matrix
PowerNetworkMatrices.DegreeTwoReduction — Type
DegreeTwoReduction <: NetworkReductionNetwork reduction algorithm that eliminates buses with exactly two connections by combining the incident branches into a single equivalent branch. This reduction preserves the electrical characteristics of the network while simplifying its topology.
Fields
irreducible_buses::Vector{Int}: List of bus numbers that should not be eliminated even if they have degree tworeduce_reactive_power_injectors::Bool: Whether to reduce buses with reactive power injections (default: true)
Examples
# Create degree-two reduction with default settings
reduction = DegreeTwoReduction()
# Create degree-two reduction protecting specific buses
reduction = DegreeTwoReduction(irreducible_buses=[101, 205])
# Create reduction that preserves buses with reactive power injections
reduction = DegreeTwoReduction(reduce_reactive_power_injectors=false)
# Apply to system
ybus = Ybus(system; network_reductions=[reduction])PowerNetworkMatrices.IncidenceMatrix — Type
Structure containing the network incidence matrix and related topology data.
The incidence matrix A represents the bus-branch connectivity of the power network, where each row corresponds to a branch and each column corresponds to a bus. Elements are:
- +1 for the "from" bus of a branch
- -1 for the "to" bus of a branch
- 0 for buses not connected to the branch
Fields
data::SparseArrays.SparseMatrixCSC{Int8, Int}: The incidence matrix data with dimensions (nbranches × nbuses). Values are {-1, 0, +1} representing the directed connectivity between branches and busesaxes::Ax: Tuple containing (arcidentifiers, busnumbers) where arcs are branch endpoint pairs and buses are the network bus numberslookup::L <: NTuple{2, Dict}: Tuple of dictionaries providing fast lookup from arc/bus identifiers to matrix indicessubnetwork_axes::Dict{Int, Ax}: Mapping from reference bus numbers to their corresponding subnetwork axesnetwork_reduction_data::NetworkReductionData: Container for network reduction information applied during matrix construction
Mathematical Properties
- Matrix Dimensions: (nbranches × nbuses)
- Element Values: {-1, 0, +1} representing directed branch-bus connectivity
- Row Sum: Each row sums to zero (conservation at branch level)
- Rank: Rank is (nbuses - nislands) for connected networks
- Sparsity: Very sparse with exactly 2 non-zero elements per branch row
Applications
- Power Flow: Forms the basis for DC power flow equations: P = A^T * f
- Sensitivity Analysis: Used in PTDF and LODF calculations
- Network Analysis: Identifies connected components and network structure
- Topology Processing: Enables network reduction and equivalencing algorithms
Notes
- Each branch contributes exactly one row with two non-zero entries (+1, -1)
- Reference buses are preserved in the matrix but identified separately
- Supports various network reduction techniques for computational efficiency
- Essential building block for BAMatrix and ABAMatrix constructions
PowerNetworkMatrices.IncidenceMatrix — Method
IncidenceMatrix(sys; network_reductions, kwargs...)
IncidenceMatrix(sys::PSY.System; network_reductions::Vector{NetworkReduction} = NetworkReduction[], kwargs...)Construct an IncidenceMatrix from a PowerSystems.System by extracting the network topology and creating the bus-branch connectivity matrix fundamental to power system analysis.
Arguments
sys::PSY.System: The power system from which to construct the incidence matrix
Keyword Arguments
network_reductions::Vector{NetworkReduction} = NetworkReduction[]: Vector of network reduction algorithms to apply before matrix constructioninclude_constant_impedance_loads::Bool=true: Whether to include constant impedance loads as shunt admittances in the network modelsubnetwork_algorithm=iterative_union_find: Algorithm used for identifying electrical islands and connected components- Additional keyword arguments are passed to the underlying
Ybusconstructor
Returns
IncidenceMatrix: The constructed incidence matrix structure containing:- Bus-branch connectivity matrix with {-1, 0, +1} elements
- Network topology information and reference bus identification
- Support for network reductions and connected component analysis
Mathematical Construction
- Network Extraction: Identifies all branches and buses from the power system
- Connectivity Mapping: Creates directed branch-bus relationships
- Matrix Assembly: Constructs sparse matrix with +1/-1 entries for branch endpoints
- Topology Analysis: Identifies reference buses and connected components
- Network Reductions: Applies specified reduction algorithms if provided
Applications
- Foundation Matrix: Essential for constructing BAMatrix and ABAMatrix
- DC Power Flow: Enables linearized power flow analysis through P = A^T * f
- Sensitivity Analysis: Required for PTDF, LODF, and other sensitivity calculations
- Network Analysis: Supports topology processing and network equivalencing
Notes
- Each branch creates exactly one matrix row with two non-zero entries
- Network reductions can significantly improve computational efficiency
- Reference buses are automatically identified for later matrix operations
- The matrix preserves full network topology for comprehensive power system analysis
PowerNetworkMatrices.IncidenceMatrix — Method
IncidenceMatrix(ybus)
IncidenceMatrix(ybus::Ybus)Construct an IncidenceMatrix from an existing Ybus matrix by extracting the network topology and creating the bus-branch connectivity matrix. This constructor leverages the network structure already captured in the Ybus matrix.
Arguments
ybus::Ybus: The Ybus matrix containing network topology and admittance data
Returns
IncidenceMatrix: The constructed incidence matrix structure containing:- Bus-branch connectivity matrix with {-1, 0, +1} elements representing directed connections
- Network topology information extracted from the Ybus structure
- Reference bus identification and subnetwork axes from the source matrix
- Network reduction data inherited from the Ybus matrix
Construction Process
- Topology Extraction: Retrieves bus and branch information from the Ybus matrix
- Arc Processing: Creates directed arc representations from branch connectivity
- Matrix Assembly: Constructs sparse incidence matrix with +1 (from bus) and -1 (to bus) entries
- Isolated Bus Handling: Includes isolated buses with zero entries for completeness
- Metadata Transfer: Preserves reference bus positions and network reduction information
Mathematical Properties
- Matrix Form: A[i,j] = +1 if branch i originates at bus j, -1 if it terminates at bus j, 0 otherwise
- Dimensions: (nbranches × nbuses) including all network branches and buses
- Sparsity: Exactly 2 non-zero entries per branch row (except for isolated buses)
- Consistency: Maintains the same network topology and reduction state as the source Ybus
Notes
- This constructor is more efficient when a Ybus matrix is already available
- Preserves all network reduction information from the source matrix
- Isolated buses are handled explicitly to maintain network completeness
- Essential for creating downstream matrices like BAMatrix and ABAMatrix from existing Ybus
PowerNetworkMatrices.LODF — Type
Structure containing the Line Outage Distribution Factor (LODF) matrix and related power system data.
The LODF matrix contains sensitivity coefficients that quantify how the outage of one transmission line affects the power flows on all other lines in the system. Each element LODF[i,j] represents the change in flow on line i when line j is taken out of service, normalized by the pre-outage flow on line j.
Fields
data::M <: AbstractArray{Float64, 2}: The LODF matrix data stored in transposed form for computational efficiency. Element (i,j) represents the sensitivity of line j flow to line i outageaxes::Ax: Tuple of identical branch/arc identifier vectors for both matrix dimensionslookup::L <: NTuple{2, Dict}: Tuple of identical dictionaries providing fast lookup from branch identifiers to matrix indicessubnetwork_axes::Dict{Int, Ax}: Mapping from reference bus numbers to their corresponding subnetwork branch axestol::Base.RefValue{Float64}: Tolerance threshold used for matrix sparsification (elements below this value are dropped)network_reduction_data::NetworkReductionData: Container for network reduction information applied during matrix construction
Mathematical Properties
- Matrix Form: LODF[i,j] = ∂fi/∂Pj where fi is flow on line i, Pj is injection change due to line j outage
- Dimensions: (nbranches × nbranches) for all transmission lines in the system
- Diagonal Elements: Always -1 (100% flow reduction on the outaged line itself)
- Symmetry: Generally non-symmetric matrix reflecting directional flow sensitivities
- Physical Meaning: Values represent fraction of pre-outage flow that redistributes to other lines
Applications
- Contingency Analysis: Evaluate impact of single line outages on system flows
- Security Assessment: Identify critical transmission bottlenecks and vulnerable lines
- System Planning: Analyze network robustness and redundancy requirements
- Real-time Operations: Support operator decision-making for preventive/corrective actions
Computational Notes
- Storage: Matrix stored in transposed form for efficient column-wise access patterns
- Sparsification: Small elements removed based on tolerance to reduce memory usage
- Linear Approximation: Based on DC power flow assumptions (neglects voltage magnitudes and reactive power)
- Single Contingencies: Designed for single line outage analysis (N-1 contingencies)
Usage Notes
- Access via
lodf[monitored_line, outaged_line]returns sensitivity coefficient - Diagonal elements are always -1.0 representing complete flow loss on outaged line
- Matrix sparsification improves performance but may introduce small numerical errors
- Results valid under DC power flow assumptions and normal operating conditions
PowerNetworkMatrices.LODF — Method
LODF(A, ABA, BA; linear_solver, tol)
LODF(A::IncidenceMatrix, ABA::ABA_Matrix, BA::BA_Matrix; linear_solver::String = "KLU", tol::Float64 = eps())Construct a Line Outage Distribution Factor (LODF) matrix from incidence, ABA, and BA matrices. This constructor provides direct control over the underlying matrix computations and is most efficient when the prerequisite matrices with factorization are already available.
Arguments
A::IncidenceMatrix: The incidence matrix containing bus-branch connectivity informationABA::ABA_Matrix: The bus susceptance matrix (A^T * B * A), preferably with KLU factorizationBA::BA_Matrix: The branch susceptance weighted incidence matrix (B * A)
Keyword Arguments
linear_solver::String = "KLU": Linear solver algorithm for matrix computations. Currently only "KLU" is supportedtol::Float64 = eps(): Sparsification tolerance for dropping small matrix elements
Returns
LODF: The constructed LODF matrix structure with line outage sensitivity coefficients
Mathematical Computation
This method computes LODF using the factorized form:
LODF = (A * ABA^(-1) * BA) / (1 - diag(A * ABA^(-1) * BA))where:
- A is the incidence matrix
- ABA^(-1) uses the factorized form from the ABA matrix (requires
ABA.Kto be factorized) - BA is the susceptance-weighted incidence matrix
Requirements and Limitations
- Factorization Required: The ABA matrix should be pre-factorized (contains KLU factorization) for efficiency
- Single Slack Bus: This method does not support distributed slack bus configurations
- Network Consistency: All three input matrices must have equivalent network reduction states
- Solver Limitation: Currently only supports "KLU" linear solver
Performance Advantages
- Pre-factorization: Leverages existing KLU factorization in ABA matrix for maximum efficiency
- Direct Computation: Avoids intermediate PTDF calculation, reducing computational steps
- Memory Efficient: Works directly with sparse matrix structures throughout computation
- Numerical Stability: Uses numerically stable KLU solver for matrix operations
Error Handling
- Validates network reduction consistency across all three input matrices
- Raises error if matrices have mismatched reduction states
- Validates linear solver selection (currently only "KLU" supported)
Usage Recommendations
- Use this constructor when you have pre-computed and factorized matrices available
- Ensure ABA matrix is factorized using
factorize(ABA)or constructed withfactorize=true - For systems with distributed slack, use the PTDF-based constructor instead
- Most efficient option for repeated LODF computations on the same network topology
PowerNetworkMatrices.LODF — Method
LODF(A, PTDFm; linear_solver, tol)
LODF(A::IncidenceMatrix, PTDFm::PTDF; linear_solver::String = "KLU", tol::Float64 = eps())Construct a Line Outage Distribution Factor (LODF) matrix from existing incidence and PTDF matrices. This constructor is more efficient when the prerequisite matrices are already available.
Arguments
A::IncidenceMatrix: The incidence matrix containing bus-branch connectivity informationPTDFm::PTDF: The power transfer distribution factor matrix (should be non-sparsified for accuracy)
Keyword Arguments
linear_solver::String = "KLU": Linear solver algorithm for matrix computations. Options: "KLU", "Dense", "MKLPardiso"tol::Float64 = eps(): Sparsification tolerance for the LODF matrix (not applied to input PTDF)
Returns
LODF: The constructed LODF matrix structure with line outage sensitivity coefficients
Mathematical Computation
The LODF matrix is computed using the formula:
LODF = (A * PTDF) / (1 - diag(A * PTDF))where:
- A is the incidence matrix representing bus-branch connectivity
- PTDF contains power transfer distribution factors
- The denominator (1 - diagonal terms) accounts for the outaged line's own flow
Important Notes
- PTDF Sparsification: The input PTDF matrix should be non-sparsified (constructed with default tolerance) to avoid accuracy issues
- Tolerance Application: The
tolparameter only affects LODF sparsification, not the input PTDF - Network Consistency: Both input matrices must have equivalent network reduction states
- Diagonal Elements: Automatically set to -1.0 representing complete flow loss on outaged lines
Performance Considerations
- Matrix Validation: Warns if input PTDF was sparsified and converts to dense format for accuracy
- Memory Usage: Sparsification with
tol > eps()can significantly reduce memory requirements - Computational Efficiency: More efficient than system-based constructor when matrices exist
Error Handling
- Validates that incidence and PTDF matrices have consistent network reduction data
- Issues warnings if sparsified PTDF matrices are used (potential accuracy issues)
- Supports automatic conversion of sparse PTDF to dense format when necessary
Linear Solver Selection
- "KLU": Recommended for most applications (sparse, numerically stable)
- "Dense": Faster for smaller systems but higher memory usage
- "MKLPardiso": Best performance for very large systems (requires MKL library)
PowerNetworkMatrices.LODF — Method
LODF(sys; linear_solver, tol, network_reductions, kwargs...)
LODF(sys::PSY.System; linear_solver::String = "KLU", tol::Float64 = eps(), network_reductions::Vector{NetworkReduction} = NetworkReduction[], kwargs...)Construct a Line Outage Distribution Factor (LODF) matrix from a PowerSystems.System by computing the sensitivity of line flows to single line outages. This is the primary constructor for LODF analysis starting from system data.
Arguments
sys::PSY.System: The power system from which to construct the LODF matrix
Keyword Arguments
linear_solver::String = "KLU": Linear solver algorithm for matrix computations. Options: "KLU", "Dense", "MKLPardiso"tol::Float64 = eps(): Sparsification tolerance for dropping small matrix elements to reduce memory usagenetwork_reductions::Vector{NetworkReduction} = NetworkReduction[]: Vector of network reduction algorithms to apply before matrix constructioninclude_constant_impedance_loads::Bool=true: Whether to include constant impedance loads as shunt admittances in the network modelsubnetwork_algorithm=iterative_union_find: Algorithm used for identifying electrical islands and connected components- Additional keyword arguments are passed to the underlying matrix constructors
Returns
LODF: The constructed LODF matrix structure containing:- Line-to-line outage sensitivity coefficients
- Network topology information and branch identifiers
- Sparsification tolerance and computational metadata
Construction Process
- Ybus Construction: Creates system admittance matrix with specified reductions
- Incidence Matrix: Builds bus-branch connectivity matrix A
- BA Matrix: Computes branch susceptance weighted incidence matrix
- PTDF Calculation: Derives power transfer distribution factors
- LODF Computation: Calculates line outage distribution factors from PTDF
- Sparsification: Applies tolerance threshold to reduce matrix density
Linear Solver Options
- "KLU": Sparse LU factorization (default, recommended for most cases)
- "Dense": Dense matrix operations (faster for small systems)
- "MKLPardiso": Intel MKL Pardiso solver (requires MKL, best for very large systems)
Mathematical Foundation
The LODF matrix is computed using the relationship:
LODF = (A * PTDF) / (1 - diag(A * PTDF))where A is the incidence matrix and PTDF is the power transfer distribution factor matrix.
Notes
- Sparsification with
tol > eps()can significantly reduce memory usage - Network reductions can improve computational efficiency for large systems
- Results are valid under DC power flow assumptions (linear approximation)
- Diagonal elements are always -1.0 representing complete flow loss on outaged lines
- For very large systems, consider using "MKLPardiso" solver with appropriate chunk size
PowerNetworkMatrices.NetworkReduction — Type
NetworkReductionAbstract base type for all network reduction algorithms used in power network analysis. Network reductions are mathematical transformations that eliminate buses and branches while preserving the electrical behavior of the remaining network elements.
Concrete implementations include:
RadialReduction: Eliminates radial (dangling) buses and branchesDegreeTwoReduction: Eliminates buses with exactly two connectionsWardReduction: Reduces external buses while preserving study bus behavior
PowerNetworkMatrices.NetworkReductionData — Type
NetworkReductionDataMutable struct containing all data mappings and metadata for network reduction operations. This structure tracks how buses and branches are mapped, combined, or eliminated during network reduction algorithms.
Fields
irreducible_buses::Set{Int}: Buses that cannot be reducedbus_reduction_map::Dict{Int, Set{Int}}: Maps retained buses to sets of eliminated busesreverse_bus_search_map::Dict{Int, Int}: Maps eliminated buses to their parent busesdirect_branch_map::Dict{Tuple{Int, Int}, PSY.ACTransmission}: One-to-one branch mappingsreverse_direct_branch_map::Dict{PSY.ACTransmission, Tuple{Int, Int}}: Reverse direct mappingsparallel_branch_map::Dict{Tuple{Int, Int}, BranchesParallel}: Parallel branch combinationsreverse_parallel_branch_map::Dict{PSY.ACTransmission, Tuple{Int, Int}}: Reverse parallel mappingsseries_branch_map::Dict{Tuple{Int, Int}, BranchesSeries}: Series branch combinationsreverse_series_branch_map::Dict{Any, Tuple{Int, Int}}: Reverse series mappingstransformer3W_map::Dict{Tuple{Int, Int}, ThreeWindingTransformerWinding}: Three-winding transformer mappingsreverse_transformer3W_map::Dict{ThreeWindingTransformerWinding, Tuple{Int, Int}}: Reverse transformer mappingsremoved_buses::Set{Int}: Set of buses eliminated from the networkremoved_arcs::Set{Tuple{Int, Int}}: Set of arcs eliminated from the networkadded_admittance_map::Dict{Int, Complex{Float32}}: Admittances added to buses during reductionadded_branch_map::Dict{Tuple{Int, Int}, Complex{Float32}}: New branches created during reductionall_branch_maps_by_type::Dict{String, Any}: Branch mappings organized by component typereductions::ReductionContainer: Container tracking applied reduction algorithmsname_to_arc_map::Dict{Type, DataStructures.SortedDict{String, Tuple{Tuple{Int, Int}, String}}}: Lazily filled with the call topopulate_branch_maps_by_type!, maps string names to their corresponding arcs and the map where the arc can be found. Used in optimization models or power flow reporting after reductions are applied. It is possible to have repeated arcs for some names if case of serial or parallel combinations.filters_applied::Dict{Type, Function}: Filters applied when populating branch maps by typedirect_branch_name_map::Dict{String, Tuple{Int, Int}}: Lazily filled, maps branch names to their corresponding arc tuples for direct branches
PowerNetworkMatrices.PTDF — Type
Structure containing the Power Transfer Distribution Factor (PTDF) matrix and related power system data.
The PTDF matrix contains sensitivity coefficients that quantify how power injections at buses affect the power flows on transmission lines. Each element PTDF[i,j] represents the incremental change in flow on line i due to a unit power injection at bus j, under DC power flow assumptions.
Fields
data::M <: AbstractArray{Float64, 2}: The PTDF matrix data stored in transposed form for computational efficiency. Element (i,j) represents the sensitivity of line j flow to bus i injectionaxes::Ax: Tuple containing (busnumbers, branchidentifiers) for matrix dimensionslookup::L <: NTuple{2, Dict}: Tuple of dictionaries providing fast lookup from bus/branch identifiers to matrix indicessubnetwork_axes::Dict{Int, Ax}: Mapping from reference bus numbers to their corresponding subnetwork axestol::Base.RefValue{Float64}: Tolerance threshold used for matrix sparsification (elements below this value are dropped)network_reduction_data::NetworkReductionData: Container for network reduction information applied during matrix construction
Mathematical Properties
- Matrix Form: PTDF[i,j] = ∂fi/∂Pj where fi is flow on line i, Pj is injection at bus j
- Dimensions: (nbuses × narcs) for all buses and impedance arcs
- Linear Superposition: Total flow = Σ(PTDF[i,j] × Pj) for all injections Pj
- Physical Meaning: Values represent the fraction of bus injection that flows through each line
- Reference Bus: Rows corresponding to reference buses are typically zero
Applications
- Power Flow Analysis: Rapid calculation of line flows for given injection patterns
- Sensitivity Studies: Evaluate impact of generation/load changes on transmission flows
- Congestion Management: Identify lines affected by specific injection changes
- Market Analysis: Support nodal pricing and transmission rights calculations
- Planning Studies: Assess transmission utilization under various scenarios
Computational Features
- Matrix Storage: Stored in transposed form (bus × branch) for efficient computation
- Sparsification: Small elements removed based on tolerance to reduce memory usage
- Reference Bus Handling: Reference bus injections automatically handled in calculations
- Distributed Slack: Supports distributed slack bus configurations for improved realism
Usage Notes
- Access via
ptdf[bus, line]returns the sensitivity coefficient - Matrix indexing uses bus numbers and branch identifiers
- Sparsification improves memory efficiency but may introduce small numerical errors
- Results valid under DC power flow assumptions (neglects voltage magnitudes and reactive power)
- Reference bus choice affects the specific values but not the relative sensitivities
PowerNetworkMatrices.PTDF — Method
PTDF(filename)
Deserialize a PTDF from an HDF5 file.
Arguments
filename::AbstractString: File containing a serialized PTDF.
PowerNetworkMatrices.PTDF — Method
PTDF(A, BA; dist_slack, linear_solver, tol)
PTDF(A::IncidenceMatrix, BA::BA_Matrix; dist_slack::Dict{Int, Float64} = Dict{Int, Float64}(), linear_solver = "KLU", tol::Float64 = eps())Construct a Power Transfer Distribution Factor (PTDF) matrix from existing incidence and BA matrices. This constructor is more efficient when the prerequisite matrices are already available and provides direct control over the underlying matrix computations.
Arguments
A::IncidenceMatrix: The incidence matrix containing bus-branch connectivity informationBA::BA_Matrix: The branch susceptance weighted incidence matrix (B × A)
Keyword Arguments
dist_slack::Dict{Int, Float64} = Dict{Int, Float64}(): Dictionary mapping bus numbers to distributed slack participation factors. Empty dictionary uses single slack bus (reference bus from matrices)linear_solver::String = "KLU": Linear solver algorithm for matrix computations. Options: "KLU", "Dense", "MKLPardiso"tol::Float64 = eps(): Sparsification tolerance for dropping small matrix elements to reduce memory usage
Returns
PTDF: The constructed PTDF matrix structure with injection-to-flow sensitivity coefficients
Mathematical Computation
The PTDF matrix is computed using the relationship:
PTDF = (A^T × B × A)^(-1) × A^T × Bwhere:
- A is the incidence matrix representing bus-branch connectivity
- B is the diagonal susceptance matrix (embedded in BA matrix)
- The computation involves solving the ABA linear system for efficiency
Distributed Slack Handling
- Single Slack: Uses reference bus identified from input matrices
- Distributed Slack: Applies participation factor corrections to final PTDF
- Automatic Processing: Dictionary converted to vector form matching matrix dimensions
- Validation: Ensures distributed slack bus numbers exist in the network
- Normalization: Participation factors automatically normalized to maintain power balance
Network Consistency Requirements
- Reduction Compatibility: Both input matrices must have equivalent network reduction states
- Reference Alignment: BA matrix reference buses determine the PTDF reference framework
- Topology Consistency: Matrices must represent the same network topology
Performance Considerations
- Matrix Reuse: More efficient when A and BA matrices are already computed
- Memory Management: Sparsification reduces storage requirements significantly
- Solver Selection: KLU recommended for sparse systems, Dense for small networks
- Computational Efficiency: Avoids redundant system matrix construction
Error Handling and Validation
- Matrix Compatibility: Validates that A and BA have consistent network reductions
- Slack Validation: Checks that distributed slack buses exist in the matrix structure
- Solver Validation: Ensures selected linear solver is supported and available
- Numerical Stability: Handles singular systems and provides informative error messages
Usage Recommendations
- Preferred Method: Use when incidence and BA matrices are already available
- Repeated Calculations: Ideal for multiple PTDF computations with different slack configurations
- Large Systems: Consider sparsification for memory efficiency
- Distributed Slack: Provides more realistic modeling of generator response to load changes
PowerNetworkMatrices.PTDF — Method
PTDF(sys; dist_slack, linear_solver, tol, kwargs...)
PTDF(sys::PSY.System; dist_slack::Dict{Int, Float64} = Dict{Int, Float64}(), linear_solver = "KLU", tol::Float64 = eps(), network_reductions::Vector{NetworkReduction} = NetworkReduction[], kwargs...)Construct a Power Transfer Distribution Factor (PTDF) matrix from a PowerSystems.System by computing the sensitivity of transmission line flows to bus power injections. This is the primary constructor for PTDF analysis starting from system data.
Arguments
sys::PSY.System: The power system from which to construct the PTDF matrix
Keyword Arguments
dist_slack::Dict{Int, Float64} = Dict{Int, Float64}(): Dictionary mapping bus numbers to distributed slack weights for realistic slack modeling. Empty dictionary uses single slack bus (default behavior)linear_solver::String = "KLU": Linear solver algorithm for matrix computations. Options: "KLU", "Dense", "MKLPardiso"tol::Float64 = eps(): Sparsification tolerance for dropping small matrix elements to reduce memory usagenetwork_reductions::Vector{NetworkReduction} = NetworkReduction[]: Vector of network reduction algorithms to apply before matrix constructioninclude_constant_impedance_loads::Bool=true: Whether to include constant impedance loads as shunt admittances in the network modelsubnetwork_algorithm=iterative_union_find: Algorithm used for identifying electrical islands and connected components- Additional keyword arguments are passed to the underlying matrix constructors
Returns
PTDF: The constructed PTDF matrix structure containing:- Bus-to-impedance-arc injection sensitivity coefficients
- Network topology information and reference bus identification
- Sparsification tolerance and computational metadata
Construction Process
- Ybus Construction: Creates system admittance matrix with specified reductions
- Incidence Matrix: Builds bus-branch connectivity matrix A
- BA Matrix: Computes branch susceptance weighted incidence matrix
- PTDF Computation: Calculates power transfer distribution factors using A^T × B^(-1) × A
- Distributed Slack: Applies distributed slack correction if specified
- Sparsification: Removes small elements based on tolerance threshold
Distributed Slack Configuration
- Single Slack: Empty
dist_slackdictionary uses conventional single slack bus - Distributed Slack: Dictionary maps bus numbers to participation factors
- Normalization: Participation factors automatically normalized to sum to 1.0
- Physical Meaning: Distributed slack better represents generator response to load changes
Linear Solver Options
- "KLU": Sparse LU factorization (default, recommended for most cases)
- "Dense": Dense matrix operations (faster for small systems, higher memory usage)
- "MKLPardiso": Intel MKL Pardiso solver (requires MKL library, best for very large systems)
Mathematical Foundation
The PTDF matrix is computed as:
PTDF = A^T × (A^T × B × A)^(-1) × A^T × Bwhere A is the incidence matrix and B is the susceptance matrix.
Notes
- Results are valid under DC power flow assumptions (linear approximation)
- Reference bus selection affects specific values but not relative sensitivities
- Sparsification with
tol > eps()can significantly reduce memory usage - Network reductions improve computational efficiency for large systems
- Distributed slack provides more realistic representation of system response
PowerNetworkMatrices.RadialReduction — Type
RadialReduction <: NetworkReductionNetwork reduction algorithm that eliminates radial (dangling) buses and their associated branches from the power network. Radial buses are leaf nodes with only one connection that do not affect the electrical behavior of the rest of the network.
Fields
irreducible_buses::Vector{Int}: List of bus numbers that should not be eliminated even if they are radial
Examples
# Create radial reduction with no protected buses
reduction = RadialReduction()
# Create radial reduction protecting specific buses
reduction = RadialReduction(irreducible_buses=[101, 205])
# Apply to system
ybus = Ybus(system; network_reductions=[reduction])PowerNetworkMatrices.VirtualLODF — Type
The Virtual Line Outage Distribution Factor (VirtualLODF) structure gathers the rows of the LODF matrix as they are evaluated on-the-go. These rows are evaluated 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.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.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.network_reduction::NetworkReduction: Structure containing the details of the network reduction applied when computing the matrix
PowerNetworkMatrices.VirtualLODF — Method
VirtualLODF(
sys;
dist_slack,
tol,
max_cache_size,
persistent_arcs,
network_reductions,
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
network_reduction::NetworkReduction: Structure containing the details of the network reduction applied when computing the matrixkwargs...: other keyword arguments used by VirtualPTDF
PowerNetworkMatrices.VirtualPTDF — Type
The Virtual Power Transfer Distribution Factor (VirtualPTDF) structure gathers the rows of the PTDF matrix as they are evaluated on-the-go. These rows are evaluated 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::Union{KLU.KLUFactorization{Float64, Int}, AppleAccelerate.AAFactorization{Float64}}: LU factorization matrices of the ABA matrix, evaluated by means of KLU or AppleAccelerateBA::SparseArrays.SparseMatrixCSC{Float64, Int}: BA matrixref_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.network_reduction::NetworkReduction: Structure containing the details of the network reduction applied when computing the matrix
PowerNetworkMatrices.VirtualPTDF — Method
VirtualPTDF(
sys;
dist_slack,
linear_solver,
tol,
max_cache_size,
persistent_arcs,
network_reductions,
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 buses.linear_solver::String = "KLU": Linear solver to use for factorization. Options: "KLU", "AppleAccelerate"tol::Float64 = eps(): Tolerance related to sparsification and values to drop.max_cache_size::Int: max cache size in MiB (initialized as MAXCACHESIZE_MiB).persistent_lines::Vector{String}: line to be evaluated as soon as the VirtualPTDF is created (initialized as empty vector of strings).network_reduction::NetworkReduction: Structure containing the details of the network reduction applied when computing the matrixkwargs...: other keyword arguments used by VirtualPTDF
PowerNetworkMatrices.Ybus — Type
Ybus{Ax, L <: NTuple{2, Dict}} <: PowerNetworkMatrix{ComplexF32}Nodal admittance matrix (Y-bus) representing the electrical admittance relationships between buses in a power system. This N×N sparse complex matrix encodes the network topology and electrical parameters needed for power flow calculations and network analysis.
Fields
data::SparseArrays.SparseMatrixCSC{ComplexF32, Int}: Sparse Y-bus matrix with complex admittance valuesadjacency_data::SparseArrays.SparseMatrixCSC{Int8, Int}: Network connectivity informationaxes::Ax: Tuple of bus axis vectors for indexing (busnumbers, busnumbers)lookup::L: Tuple of lookup dictionaries mapping bus numbers to matrix indicessubnetwork_axes::Dict{Int, Ax}: Bus axes for each electrical island/subnetworkarc_subnetwork_axis::Dict{Int, Vector{Tuple{Int, Int}}}: Arc axes for each subnetworknetwork_reduction_data::NetworkReductionData: Metadata from network reduction operationsarc_admittance_from_to::Union{ArcAdmittanceMatrix, Nothing}: From-to arc admittance matrixarc_admittance_to_from::Union{ArcAdmittanceMatrix, Nothing}: To-from arc admittance matrix
Key Features
- Indexed by bus numbers (non-sequential numbering supported)
- Supports network reductions (radial, degree-two, Ward)
- Handles multiple electrical islands/subnetworks
- Optional arc admittance matrices for power flow calculations
- Sparse matrix representation for computational efficiency
Usage
The Y-bus is fundamental for:
- Power flow analysis: V = Y⁻¹I
- Short circuit calculations
- Network impedance analysis
- Sensitivity analysis (PTDF/LODF)
Examples
# Basic Y-bus construction
ybus = Ybus(system)
# With arc admittance matrices for power flow
ybus = Ybus(system; make_arc_admittance_matrices=true)
# With network reductions
ybus = Ybus(system; network_reductions=[RadialReduction(), DegreeTwoReduction()])See Also
PTDF: Power Transfer Distribution FactorsLODF: Line Outage Distribution FactorsNetworkReduction: Network reduction algorithms
PowerNetworkMatrices.Ybus — Method
Ybus(
sys;
make_arc_admittance_matrices,
network_reductions,
include_constant_impedance_loads,
subnetwork_algorithm
)
Ybus(
sys::PSY.System;
make_arc_admittance_matrices::Bool = false,
network_reductions::Vector{NetworkReduction} = NetworkReduction[],
include_constant_impedance_loads::Bool = true,
subnetwork_algorithm = iterative_union_find,
kwargs...
) -> YbusConstruct a nodal admittance matrix (Y-bus) from a power system.
Builds the sparse complex Y-bus matrix representing the electrical admittance relationships between buses in the power system. Handles AC branches, transformers, shunt elements, and network reductions while maintaining connectivity analysis.
Arguments
sys::PSY.System: Power system to build Y-bus from
Keyword arguments
make_arc_admittance_matrices::Bool=false: Whether to construct arc admittance matrices for power flownetwork_reductions::Vector{NetworkReduction}=[]: Network reduction algorithms to applyinclude_constant_impedance_loads::Bool=true: Whether to include constant impedance loads as shunt admittancessubnetwork_algorithm=iterative_union_find: Algorithm for finding electrical islands
Returns
Ybus: Constructed Y-bus matrix with network topology and electrical parameters
Features
- Branch Support: Lines, transformers, phase shifters, three-winding transformers
- Shunt Elements: Fixed admittances, switched admittances, constant impedance loads
- Network Reductions: Radial, degree-two, Ward reductions for computational efficiency
- Multiple Islands: Handles disconnected network components with separate reference buses
- Branch Matrices: Optional from-to/to-from admittance matrices for power flow calculations
Examples
# Basic Y-bus construction
ybus = Ybus(system)
# With arc admittance matrices for power flow
ybus = Ybus(system; make_arc_admittance_matrices=true)
# Apply network reductions for computational efficiency
reductions = [RadialReduction(), DegreeTwoReduction()]
ybus = Ybus(system; network_reductions=reductions)
# Exclude constant impedance loads
ybus = Ybus(system; include_constant_impedance_loads=false)See Also
NetworkReduction: Network reduction algorithmsPTDF: Power transfer distribution factorsLODF: Line outage distribution factors
PowerNetworkMatrices.depth_first_search — Method
depth_first_search(M, bus_numbers)
depth_first_search(M::SparseArrays.SparseMatrixCSC, bus_numbers::Vector{Int})Find connected subnetworks using depth-first search algorithm.
Arguments
M::SparseArrays.SparseMatrixCSC: Sparse matrix representing network connectivitybus_numbers::Vector{Int}: Vector containing the bus numbers of the system
Returns
Dict{Int, Set{Int}}: Dictionary mapping representative bus numbers to sets of connected buses
PowerNetworkMatrices.factorize — Method
factorize(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 — Method
find_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 — Method
find_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 — Method
find_subnetworks(M, bus_numbers; subnetwork_algorithm)
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.subnetwork_algorithm::Function: algorithm for computing subnetworks. Valid options are iterativeunionfind (default) and depthfirstsearch
PowerNetworkMatrices.find_subnetworks — Method
find_subnetworks(M)
find_subnetworks(M::Ybus) -> Dict{Int, Set{Int}}Identify electrical islands (subnetworks) in the Y-bus matrix.
Analyzes the network topology to find groups of buses that are electrically connected to each other but isolated from other groups. Each subnetwork represents an electrical island that requires its own reference bus and can be solved independently.
Arguments
M::Ybus: Y-bus matrix to analyze
Returns
Dict{Int, Set{Int}}: Dictionary mapping reference bus numbers to sets of bus numbers in each subnetwork
Examples
ybus = Ybus(system)
subnetworks = find_subnetworks(ybus)
for (ref_bus, buses) in subnetworks
println("Island ", ref_bus, ": ", sort(collect(buses)))
end
if length(subnetworks) > 1
@warn "Network has ", length(subnetworks), " electrical islands"
endImplementation Details
- Uses adjacency matrix analysis to find connected components
- Each subnetwork gets assigned a reference bus for voltage angle reference
- Isolated buses or groups require separate power flow analysis
- Critical for power flow initialization and solution
See Also
validate_connectivity: Check for full connectivitydepth_first_search: Graph traversal algorithmiterative_union_find: Alternative connectivity algorithm
PowerNetworkMatrices.from_hdf5 — Method
from_hdf5(_, filename)
Deserialize a PTDF from an HDF5 file.
Arguments
::Type{PTDF}:filename::AbstractString: File containing a serialized PTDF.
PowerNetworkMatrices.get_bus_reduction_map — Method
get_bus_reduction_map(rb)
get_bus_reduction_map(rb::NetworkReductionData)Get the bus reduction map from NetworkReductionData.
Arguments
rb::NetworkReductionData: The network reduction data
Returns
Dict{Int, Set{Int}}: Dictionary mapping retained buses to sets of removed buses
PowerNetworkMatrices.get_lodf_data — Method
get_lodf_data(lodf)
get_lodf_data(lodf::LODF)Extract the LODF matrix data in the standard orientation (non-transposed).
Arguments
lodf::LODF: The LODF structure from which to extract data
Returns
AbstractArray{Float64, 2}: The LODF matrix data with standard orientation
PowerNetworkMatrices.get_ptdf_data — Method
get_ptdf_data(ptdf)
get_ptdf_data(ptdf::PTDF)Extract the PTDF matrix data in the standard orientation (non-transposed).
Arguments
ptdf::PTDF: The PTDF structure from which to extract data
Returns
AbstractArray{Float64, 2}: The PTDF matrix data with standard orientation
PowerNetworkMatrices.get_reductions — Method
get_reductions(rb)
get_reductions(rb::NetworkReductionData)Get the reduction container from NetworkReductionData.
Arguments
rb::NetworkReductionData: The network reduction data
Returns
ReductionContainer: Container with the applied network reductions
PowerNetworkMatrices.get_ward_reduction — Method
get_ward_reduction(
data,
bus_lookup,
bus_axis,
boundary_buses,
ref_bus_numbers,
study_buses
)
get_ward_reduction(data, bus_lookup, bus_axis, boundary_buses, ref_bus_numbers, study_buses)Perform Ward reduction to create an equivalent network representation.
Ward reduction is a network reduction technique that eliminates external buses while preserving the electrical characteristics seen from the study buses. External buses are mapped to boundary buses based on impedance criteria, and equivalent admittances are computed.
Arguments
data::SparseArrays.SparseMatrixCSC{ComplexF32, Int}: Admittance matrix of the systembus_lookup::Dict{Int, Int}: Dictionary mapping bus numbers to matrix indicesbus_axis::Vector{Int}: Vector of all bus numbers in the systemboundary_buses::Set{Int}: Set of boundary bus numbers between study and external areasref_bus_numbers::Set{Int}: Set of reference bus numbersstudy_buses::Vector{Int}: Vector of study bus numbers to retain
Returns
Tuple: Contains bus reduction map, reverse bus search map, added branch map, and added admittance map
PowerNetworkMatrices.is_factorized — Method
is_factorized(ABA)
is_factorized(ABA::ABA_Matrix)Check if an ABA_Matrix has been factorized (i.e., contains LU factorization matrices).
Arguments
ABA::ABA_Matrix: The ABA matrix to check
Returns
Bool: true if the matrix has been factorized, false otherwise
PowerNetworkMatrices.iterative_union_find — Method
iterative_union_find(M, bus_numbers)
iterative_union_find(M::SparseArrays.SparseMatrixCSC, bus_numbers::Vector{Int})Find connected subnetworks using iterative union-find algorithm.
Arguments
M::SparseArrays.SparseMatrixCSC: Sparse matrix representing network connectivitybus_numbers::Vector{Int}: Vector containing the bus numbers of the system
Returns
Dict{Int, Set{Int}}: Dictionary mapping representative bus numbers to sets of connected buses
PowerNetworkMatrices.to_hdf5 — Method
to_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 — Method
validate_connectivity(M)
Validates connectivity by checking that the number of subnetworks is 1 (fully connected network).
PowerNetworkMatrices.validate_connectivity — Method
validate_connectivity(sys)
Checks the network connectivity of the system using Depth First Search (DFS)
PowerNetworkMatrices.validate_connectivity — Method
validate_connectivity(M)
validate_connectivity(M::Ybus) -> BoolValidate that the Y-bus represents a fully connected electrical network.
Checks network connectivity by counting the number of electrical islands (subnetworks) in the Y-bus matrix. A fully connected network should have exactly one subnetwork. Multiple subnetworks indicate electrical isolation between parts of the system.
Arguments
M::Ybus: Y-bus matrix to validate
Returns
Bool:trueif network is fully connected (single subnetwork),falseotherwise
Examples
ybus = Ybus(system)
if validate_connectivity(ybus)
println("Network is fully connected")
else
println("Network has isolated islands")
islands = find_subnetworks(ybus)
println("Number of islands: ", length(islands))
endImplementation Details
- Uses
find_subnetworks()to identify electrical islands - Single subnetwork indicates full electrical connectivity
- Multiple subnetworks may require separate power flow solutions
See Also
find_subnetworks: Identify electrical islandsvalidate_connectivity: System-level connectivity validation