Internal - Core
Power Flow Types
PowerFlowData
Struct and Type Definitions
PowerFlows.ABAPowerFlowData — Type
A type alias for a PowerFlowData struct whose type parameters are configured for the DCPowerFlow method.
PowerFlows.ACPowerFlowData — Type
A type alias for a PowerFlowData struct whose type parameters are configured for the ACPowerFlow method.
PowerFlows.PTDFPowerFlowData — Type
A type alias for a PowerFlowData struct whose type parameters are configured for the PTDFDCPowerFlow method .
PowerFlows.PowerFlowData — Type
PowerFlowData{M <: PNM.PowerNetworkMatrix, N <: Union{PNM.PowerNetworkMatrix, Nothing}}Structure containing all the data required for the evaluation of the power flows and angles, as well as these ones.
All fields starting with bus_ are ordered according to bus_lookup, and all fields starting with arc_ are ordered according to arc_lookup: one row per bus/arc, one column per time period. Here, buses should be understood as "buses remaining, after the network reduction." Similarly, we use "arcs" instead of "branches" to distinguish between network elements (post-reduction) and system objects (pre-reduction).
Generally, do not construct this directly. Instead, use one of the later constructors to pass in a PowerFlowEvaluationModel and a PowerSystems.System. aux\_network\_matrix and power\_network\_matrix will then be set to the appropriate matrices that are needed for computing that type of power flow. See also ACPowerFlowData, ABAPowerFlowData, PTDFPowerFlowData, and vPTDFPowerFlowData: these are all aliases for PowerFlowData{N, M} with specific N,M, that are used for the respective type of power flow evaluations.
Fields:
bus_active_power_injections::Matrix{Float64}: matrix containing the bus active power injections.bus_reactive_power_injections::Matrix{Float64}: matrix containing the bus reactive power injections.bus_active_power_withdrawals::Matrix{Float64}: matrix containing the bus reactive power withdrawals.bus_reactive_power_withdrawals::Matrix{Float64}: matrix containing the bus reactive power withdrawals.bus_active_power_constant_current_withdrawals::Matrix{Float64}: matrix containing the bus active power constant current withdrawals.bus_reactive_power_constant_current_withdrawals::Matrix{Float64}: matrix containing the bus reactive power constant current withdrawals.bus_active_power_constant_impedance_withdrawals::Matrix{Float64}: matrix containing the bus active power constant impedance withdrawals.bus_reactive_power_constant_impedance_withdrawals::Matrix{Float64}: matrix containing the bus reactive power constant impedance withdrawals.bus_reactive_power_bounds::Matrix{Float64}: matrix containing upper and lower bounds for the reactive supply at each bus at each time period.bus_type::Matrix{PSY.ACBusTypes}: matrix containing type of buses present in the system.bus_magnitude::Matrix{Float64}: matrix containing the bus voltage magnitudes.bus_angles::Matrix{Float64}: matrix containing the bus voltage angles.arc_active_power_flow_from_to::Matrix{Float64}: matrix containing the active power flows measured at thefrombus.arc_reactive_power_flow_from_to::Matrix{Float64}: matrix containing the reactive power flows measured at thefrombus.arc_active_power_flow_to_from::Matrix{Float64}: matrix containing the active power flows measured at thetobus.arc_reactive_power_flow_to_from::Matrix{Float64}: matrix containing the reactive power flows measured at thetobus.arc_angle_differences::Matrix{Float64}: matrix containing the voltage angle difference (θfrom − θto) across each arc.generic_hvdc_flows::Dict{Tuple{Int, Int}, Tuple{Float64, Float64}}: dictionary mapping each generic HVDC line (represented as a tuple of the from and to bus numbers) to a tuple of(P_from_to, P_to_from)active power flows.bus_hvdc_net_power::Matrix{Float64}: "(b, t)" matrix containing the net power injections from all HVDC lines at each bus. b: number of buses, t: number of time period. Only contains HVDCs handled as separate injection/withdrawal pairs: LCCs and generic for DC, or just generic for AC.time_step_map::Dict{Int, S}: dictionary mapping the number of the time periods (corresponding to the column number of the previously mentioned matrices) and their names.power_network_matrix::M: matrix used for the evaluation of either the power flows or bus angles, depending on the method considered.aux_network_matrix::N: matrix used for the evaluation of either the power flows or bus angles, depending on the method considered.neighbors::Vector{Set{Int}}: Vector with the sets of adjacent buses.
PowerFlows.PowerFlowData — Method
PowerFlowData(
pf::ACPowerFlow{<:ACPowerFlowSolverType},
sys::PSY.System
) -> ACPowerFlowData{<:ACPowerFlowSolverType}Creates the structure for an AC power flow calculation, given the System sys. Configuration options like time_steps, timestep_names, network_reductions, and correct_bustypes are taken from the ACPowerFlow object.
Calling this function will not evaluate the power flows and angles. Note that first input is of type ACPowerFlow: this version is used to solve AC power flows, and returns an ACPowerFlowData object.
Arguments:
pf::ACPowerFlow: the settings for the AC power flow solver, includingtime_steps,time_step_names,network_reductions, andcorrect_bustypes.sys::PSY.System: ASystemobject that represents the power grid under consideration.
WARNING: functions for the evaluation of the multi-period AC PF still to be implemented.
PowerFlows.PowerFlowData — Method
PowerFlowData(
pf::DCPowerFlow,
sys::PSY.System
) -> ABAPowerFlowDataCreates a PowerFlowData structure configured for a standard DC power flow calculation, given the System sys. Configuration options like time_steps, time_step_names, network_reductions, and correct_bustypes are taken from the DCPowerFlow object.
Calling this function will not evaluate the power flows and angles. Note that first input is of type DCPowerFlow: this version is used to solve DC power flows, and returns an ABAPowerFlowData object.
Arguments:
pf::DCPowerFlow: Run a DC power flow: internally, store the ABA matrix aspower_network_matrixand the BA matrix asaux_network_matrix. Configuration options are taken from this object.sys::PSY.System: ASystemobject that represents the power grid under consideration.
PowerFlows.PowerFlowData — Method
PowerFlowData(
pf::PTDFDCPowerFlow,
sys::PSY.System
) -> PTDFPowerFlowDataCreates a PowerFlowData structure configured for a Partial Transfer Distribution Factor Matrix DC power flow calculation, given the System sys. Configuration options like time_steps, time_step_names, network_reductions, and correct_bustypes are taken from the PTDFDCPowerFlow object.
Calling this function will not evaluate the power flows and angles. Note that first input is of type PTDFDCPowerFlow: this version is used to solve DC power flows via the Power Transfer Distribution Factor (PTDF) matrix. This function returns a PTDFPowerFlowData object.
Arguments:
pf::PTDFDCPowerFlow: Run a DC power flow with PTDF matrix: internally, store the PTDF matrix aspower_network_matrixand the ABA matrix asaux_network_matrix. Configuration options are taken from this object.sys::PSY.System: ASystemobject that represents the power grid under consideration.
PowerFlows.PowerFlowData — Method
PowerFlowData(
pf::vPTDFDCPowerFlow,
sys::PSY.System
) -> vPTDFPowerFlowDataCreates a PowerFlowData structure configured for a virtual Partial Transfer Distribution Factor Matrix DC power flow calculation, given the System sys. Configuration options like time_steps, time_step_names, network_reductions, and correct_bustypes are taken from the vPTDFDCPowerFlow object.
Calling this function will not evaluate the power flows and angles. Note that first input is of type vPTDFDCPowerFlow: this version is used to solve DC power flows using a virtual Power Transfer Distribution Factor (PTDF) matrix. This function returns a vPTDFPowerFlowData object.
Arguments:
pf::vPTDFDCPowerFlow: Run a virtual PTDF power flow: internally, store the virtual PTDF matrix aspower_network_matrixand the ABA matrix asaux_network_matrix. Configuration options are taken from this object.sys::PSY.System: ASystemobject that represents the power grid under consideration.
PowerFlows.PowerFlowData — Method
Sets the two PowerNetworkMatrix fields and a few others (time_steps, time_step_map), then creates arrays of default values (usually zeros) for the rest.
PowerFlows.SystemPowerFlowContainer — Type
A PowerFlowContainer that represents its data as a PSY.System.
PowerFlows.vPTDFPowerFlowData — Type
A type alias for a PowerFlowData struct whose type parameters are configured for the vPTDFDCPowerFlow method.
PowerFlows._compute_arc_angle_differences_from_data! — Method
Compute arc angle differences for all arcs and all time steps, looking up bus indices from the arc axis and bus lookup stored in data. Used by DC solvers.
PowerFlows._compute_arc_angle_differences_from_indices! — Method
Compute arc angle differences using precomputed from/to bus index vectors over specified time steps. Used by the AC solver where fb_ix/tb_ix are already available from the branch flow calculation.
PowerFlows.make_power_flow_container — Function
Create an appropriate PowerFlowContainer for the given PowerFlowEvaluationModel and initialize it from the given PSY.System.
Configuration options like time_steps, time_step_names, network_reductions, and correct_bustypes are taken from the PowerFlowEvaluationModel object.
Arguments:
pfem::PowerFlowEvaluationModel: power flow model to construct a container for (e.g.,DCPowerFlow())sys::PSY.System: the System from which to initialize the power flow container
PowerFlows.supports_multi_period — Method
Trait signifying whether the PowerFlowContainer can represent multi-period data. Must be implemented for all concrete subtypes.
Solving a PowerFlowData instance
PowerFlows.solve_power_flow! — Method
solve_power_flow!(data::ACPowerFlowData; kwargs...)Solve the multiperiod AC power flow problem for the given power flow data.
The bus types can be changed from PV to PQ if the reactive power limits are violated. The power flow solver settings are taken from the ACPowerFlow object stored in data.
Arguments
data::ACPowerFlowData: The power flow data containing the grid information and initial conditions.kwargs...: Additional keyword arguments. If these overlap with those in thesolver_settingsof theACPowerFlowobject, the values inkwargstake precedence.
Keyword Arguments
time_steps: Specifies the time steps to solve. Defaults to sorting and collecting the keys ofget_time_step_map(data).
Description
This function solves the AC power flow problem for each time step specified in data. It preallocates memory for the results and iterates over the sorted time steps. For each time step, it calls the _ac_power_flow function to solve the power flow equations and updates the data object with the results. If the power flow converges, it updates the active and reactive power injections, as well as the voltage magnitudes and angles for different bus types (REF, PV, PQ). If the power flow does not converge, it sets the corresponding entries in data to NaN. Finally, it calculates the branch power flows and updates the data object.
Notes
- If the grid topology changes (e.g., tap positions of transformers or in-service status of branches), the admittance matrices
YftandYtfmust be updated. - If
YftandYtfchange between time steps, the branch flow calculations must be moved inside the loop.
Examples
solve_power_flow!(data)PowerFlows._get_arc_resistances — Method
_get_arc_resistances(data::Union{PTDFPowerFlowData, vPTDFPowerFlowData, ABAPowerFlowData}) -> Vector{Float64}Look up the equivalent resistance of each arc from the network reduction data, using PNM API where available.
PowerFlows.adjust_power_injections_for_lccs! — Method
Adjust the power injections vector to account for the power flows through LCCs.
Relies on the fact that we calculate those flows during initialization and save them to the active_power_flow_from_to and active_power_flow_to_from fields of the LCCParameters struct.
PowerFlows.dc_loss_factors — Method
dc_loss_factors(
data::Union{PTDFPowerFlowData, vPTDFPowerFlowData},
P::Matrix{Float64},
) -> Matrix{Float64}Compute the gradient of total system active power losses with respect to bus injections using the DC power flow approximation:
∂Loss/∂P = 2 · PTDFᵀ · diag(R) · PTDF · PThis is equivalent to the per-element form:
∂Loss/∂Pᵢ = Σₖ 2·Rₖ·PTDFₖᵢ·Σⱼ PTDFₖⱼ·PⱼArguments
data::Union{PTDFPowerFlowData, vPTDFPowerFlowData}: solved power flow data containing the PTDF matrix and network reduction data for looking up branch resistances.P::Matrix{Float64}: bus injection matrix of size(num_buses, num_timesteps).
Returns
Matrix{Float64}: loss factor matrix of size(num_buses, num_timesteps), where each entry[i, t]is the marginal change in total system losses per unit injection at busiin time stept.
PowerFlows.solve_power_flow! — Method
solve_power_flow!(data::PTDFPowerFlowData)Evaluates the PTDF power flow and writes the result to the fields of the PTDFPowerFlowData structure.
This function modifies the following fields of data, setting them to the computed values:
data.bus_angles: the bus angles for each bus in the system.data.branch_active_power_flow_from_to: the active power flow from the "from" bus to the "to" bus of each branchdata.branch_active_power_flow_to_from: the active power flow from the "to" bus to the "from" bus of each branch
Additionally, it sets data.converged to true, indicating that the power flow calculation was successful.
PowerFlows.solve_power_flow! — Method
solve_power_flow!(data::vPTDFPowerFlowData)Evaluates the virtual PTDF power flow and writes the results to the fields of the vPTDFPowerFlowData structure.
This function modifies the following fields of data, setting them to the computed values:
data.bus_angles: the bus angles for each bus in the system.data.branch_active_power_flow_from_to: the active power flow from the "from" bus to the "to" bus of each branchdata.branch_active_power_flow_to_from: the active power flow from the "to" bus to the "from" bus of each branch
Additionally, it sets data.converged to true, indicating that the power flow calculation was successful.
Manipulating a PowerFlowData instance
PowerFlows.partition_state — Method
Partitions the state vector's variables based on what physical quantity each represents. Returns a NamedTuple, with the 4 keys Va, Vm, P, and Q. The 4 values are vectors of length equal to the number of buses, with NaNs in the positions where that physical quantity is not part of the state vector for that bus. (Currently not intended for use in spots where performance is critical.)
PowerFlows.update_data! — Method
Update the fields of data based on the values of the state vector.
PowerFlows.update_net_power! — Method
Update P_net and Q_net based on the values of the state vector.
PowerFlows.update_state! — Method
Update state vector based on values of fields of data.
PowerFlows.initialize_power_flow_data! — Method
Sets the fields of a PowerFlowData struct to match the given System.
PowerFlows._dc_power_flow_fallback! — Method
When solving AC power flows, if the initial guess has large residual, we run a DC power flow as a fallback. This runs a DC power flow on data::ACPowerFlowData for the given time_step, and writes the solution to data.bus_angles.
PowerFlows.calculate_x0 — Method
Calculate x0 from data.
PowerFlows.dc_power_flow_start! — Method
If initial residual is large, run a DC power flow and see if that gives a better starting point for angles. If so, then overwrite x0 with the result of the DC power flow. If not, keep the original x0.
LCC HVDC Parameters and Utilities
PowerFlows._calculate_dQ_dV_lcc — Method
_calculate_dQ_dV_lcc(t::Float64, I_dc::Float64, x_t::Float64, Vm::Float64, ϕ::Float64) -> Float64Compute the derivative of reactive power Q with respect to voltage magnitude Vm for LCC converter calculations.
PowerFlows._calculate_dQ_dt_lcc — Method
_calculate_dQ_dt_lcc(t::Float64, I_dc::Float64, x_t::Float64, Vm::Float64, ϕ::Float64) -> Float64Compute the derivative of reactive power Q with respect to transformer tap t for LCC converter calculations.
PowerFlows._calculate_dQ_dα_lcc — Method
_calculate_dQ_dα_lcc(t::Float64, I_dc::Float64, x_t::Float64, Vm::Float64, ϕ::Float64, α::Float64) -> Float64Compute the derivative of reactive power Q with respect to firing/extinction angle α for LCC converter calculations.
PowerFlows._calculate_y_lcc — Method
_calculate_y_lcc(t::Float64, I_dc::Float64, Vm::Float64, ϕ::Float64) -> ComplexF64Compute the admittance value Y for LCC converter calculations.
PowerFlows._calculate_ϕ_lcc — Method
_calculate_ϕ_lcc(α::Float64, I_dc::Float64, x_t::Float64, Vm::Float64) -> Float64Compute the phase angle ϕ for LCC converter calculations.
PowerFlows.hvdc_fixed_injections! — Method
Adjust the power injections/withdrawal vectors to account for all HVDC lines of a given type, modeling those HVDC lines as a simple fixed injection/withdrawal at each terminal.
PowerFlows.initialize_LCC_arcs_and_buses! — Method
Initialize the arcs and bus_indices fields of the LCCParameters structure in the PowerFlowData.
AC Power Flow
Residuals
PowerFlows.ACPowerFlowResidual — Type
struct ACPowerFlowResidualA struct to keep track of the residuals in the Newton-Raphson AC power flow calculation.
Fields
data::ACPowerFlowData: The grid model data.Rf!::Function: A function that updates the residuals based on the latest values stored in the grid at the given iteration.Rv::Vector{Float64}: A vector of the values of the residuals.P_net::Vector{Float64}: A vector of net active power injections.Q_net::Vector{Float64}: A vector of net reactive power injections.P_net_set::Vector{Float64}: A vector of the set-points for active power injections (their initial values before power flow calculation).bus_slack_participation_factors::SparseVector{Float64, Int}: A sparse vector of the slack participation factors aggregated at the bus level.subnetworks::Dict{Int64, Vector{Int64}}: The dictionary that identifies subnetworks (connected components), with the key defining the REF bus, values defining the corresponding buses in the subnetwork.
PowerFlows.ACPowerFlowResidual — Method
ACPowerFlowResidual(data::ACPowerFlowData, time_step::Int64)Create an instance of ACPowerFlowResidual for a given time step.
Arguments
data::ACPowerFlowData: The power flow data representing the power system model.time_step::Int64: The time step for which the power flow calculation is executed.
Returns
ACPowerFlowResidual: An instance containing the residual values, net bus active power injections, and net bus reactive power injections.
PowerFlows.ACPowerFlowResidual — Method
(Residual::ACPowerFlowResidual)(x::Vector{Float64}, time_step::Int64)Update the AC power flow residuals inplace and store the result in the attribute Rv of the struct. The inputs are the values of state vector x and the current time step time_step. This function implements the functor approach for the ACPowerFlowResidual struct. This makes the struct callable. Calling the ACPowerFlowResidual will also update the values of P, Q, V, Θ in the data struct.
Arguments
x::Vector{Float64}: The state vector values.time_step::Int64: The current time step.
PowerFlows.ACPowerFlowResidual — Method
(Residual::ACPowerFlowResidual)(Rv::Vector{Float64}, x::Vector{Float64}, time_step::Int64)Evaluate the AC power flow residuals and store the result in Rv using the provided state vector x and the current time step time_step. The residuals are updated inplace in the struct and additionally copied to the provided array. This function implements the functor approach for the ACPowerFlowResidual struct. This makes the struct callable. Calling the ACPowerFlowResidual will also update the values of P, Q, V, Θ in the data struct.
Arguments
Rv::Vector{Float64}: The vector to store the calculated residuals.x::Vector{Float64}: The state vector.time_step::Int64: The current time step.
PowerFlows._update_residual_values! — Method
_update_residual_values!(
F::Vector{Float64},
x::Vector{Float64},
P_net::Vector{Float64},
Q_net::Vector{Float64},
data::ACPowerFlowData,
time_step::Int64,
)Update the residual values for the Newton-Raphson AC power flow calculation. This function is used internally in the ACPowerFlowResidual struct. This function also updates the values of P, Q, V, Θ in the data struct.
Arguments
F::Vector{Float64}: Vector of the values of the residuals.x::Vector{Float64}: State vector values.P_net::Vector{Float64}: Vector of net active power injections at each bus.Q_net::Vector{Float64}: Vector of net reactive power injections at each bus.P_net_set::Vector{Float64}: Vector of the set-points for active power injections (their initial values before power flow calculation).bus_slack_participation_factors::SparseVector{Float64, Int}: Sparse vector of the slack participation factors aggregated at the bus level.ref_bus::Int: The index of the reference bus to be used for the total slack power.data::ACPowerFlowData: Data structure representing the grid model for the AC power flow calculation.time_step::Int64: The current time step for which the residual values are being updated.
Jacobian
PowerFlows.ACPowerFlowJacobian — Type
struct ACPowerFlowJacobianA struct that represents the Jacobian matrix for AC power flow calculations.
This struct uses the functor pattern, meaning instances of ACPowerFlowJacobian store the data (Jacobian matrix) internally and can be called as a function at the same time. Calling the instance as a function updates the stored Jacobian matrix.
Fields
data::ACPowerFlowData: The grid model data used for power flow calculations.Jf!::Function: A function that calculates the Jacobian matrix inplace.Jv::SparseArrays.SparseMatrixCSC{Float64, Int32}: The Jacobian matrix, which is updated by the functionJf!.diag_elements::MVector{4, Float64}: Temporary storage for diagonal elements during Jacobian update.bus_slack_participation_factors::SparseVector{Float64, Int}: Normalized per-bus slack participation factors for the current time step (from theACPowerFlowResidual). Used for the distributed slack Jacobian entries.subnetworks::Dict{Int64, Vector{Int64}}: Subnetwork mapping from REF bus to bus list (from theACPowerFlowResidual). Used for the distributed slack Jacobian entries.
PowerFlows.ACPowerFlowJacobian — Method
(J::ACPowerFlowJacobian)(time_step::Int64)Update the Jacobian matrix Jv using the function Jf! and the provided data and time step.
Defining this method allows an instance of ACPowerFlowJacobian to be called as a function, following the functor pattern.
Arguments
time_step::Int64: The time step for the calculations.
Example
residual = ACPowerFlowResidual(data, time_step)
J = ACPowerFlowJacobian(data,
residual.bus_slack_participation_factors,
residual.subnetworks,
time_step
)
J(time_step) # Updates the Jacobian matrix JvPowerFlows.ACPowerFlowJacobian — Method
ACPowerFlowJacobian(data::ACPowerFlowData,
bus_slack_participation_factors::SparseVector{Float64, Int},
subnetworks::Dict{Int64, Vector{Int64}}, time_step::Int64) -> ACPowerFlowJacobianThis is the constructor for ACPowerFlowJacobian. Create an ACPowerFlowJacobian instance. As soon as the instance is created, it already has the Jacobian matrix structure initialized and its values updated, stored internally as Jv. The data instance is stored internally and used to update the Jacobian matrix because the structure of the Jacobian matrix is tied to the data. Changing the data requires creating a new instance of ACPowerFlowJacobian.
The bus_slack_participation_factors and subnetworks are needed for the distributed slack Jacobian entries: each participating bus k has ∂FPk/∂x[2*ref-1] = -c_k in the Jacobian.
Arguments
data::ACPowerFlowData: The data used for power flow calculations.bus_slack_participation_factors::SparseVector{Float64, Int}: Normalized per-bus slack participation factors (from theACPowerFlowResidual).subnetworks::Dict{Int64, Vector{Int64}}: Subnetwork mapping from REF bus to bus list (from theACPowerFlowResidual).time_step::Int64: The time step for the calculations.
Returns
ACPowerFlowJacobian: An instance ofACPowerFlowJacobian.
Example
residual = ACPowerFlowResidual(data, time_step)
J = ACPowerFlowJacobian(data,
residual.bus_slack_participation_factors, residual.subnetworks, time_step)
J(time_step) # Updates the Jacobian matrix stored internally in J.
J.Jv # Access the Jacobian matrix stored internally in J.PowerFlows.ACPowerFlowJacobian — Method
(J::ACPowerFlowJacobian)(J::SparseArrays.SparseMatrixCSC{Float64, Int32}, time_step::Int64)Use the ACPowerFlowJacobian to update the provided Jacobian matrix J inplace.
Update the internally stored Jacobian matrix Jv using the function Jf! and the provided data and time step, and write the updated Jacobian values to J.
This method allows an instance of ACPowerFlowJacobian to be called as a function, following the functor pattern.
Arguments
J::SparseArrays.SparseMatrixCSC{Float64, Int32}: A sparse matrix to be updated with new values of the Jacobian matrix.time_step::Int64: The time step for the calculations.
Example
residual = ACPowerFlowResidual(data, time_step)
J = ACPowerFlowJacobian(data,
residual.bus_slack_participation_factors,
residual.subnetworks,
time_step
)
Jv = SparseArrays.sparse(Float64[], J_INDEX_TYPE[], J_INDEX_TYPE[])
J(Jv, time_step) # Updates the Jacobian matrix Jv and writes it to JPowerFlows._block_J_indices — Method
_block_J_indices(data::ACPowerFlowData, time_step::Int) -> (Vector{Int32}, Vector{Int32})Get the indices to reindex the Jacobian matrix from the interleaved form to the block form:
\[\begin{bmatrix} \frac{\partial P}{\partial \theta} & \frac{\partial P}{\partial V} \\ \frac{\partial Q}{\partial \theta} & \frac{\partial Q}{\partial V} \end{bmatrix}\]
Arguments
pvpq::Vector{Int32}: Indices of the buses that are PV or PQ buses.pq::Vector{Int32}: Indices of the buses that are PQ buses.
Returns
rows::Vector{Int32}: Row indices for the block Jacobian matrix.cols::Vector{Int32}: Column indices for the block Jacobian matrix.
PowerFlows._calculate_loss_factors — Method
calculate_loss_factors(data::ACPowerFlowData, Jv::SparseMatrixCSC{Float64, Int32}, time_step::Int)Calculate and store the active power loss factors in the loss_factors matrix of the ACPowerFlowData structure for a given time step.
The loss factors are computed using the Jacobian matrix Jv and the vector dSbus_dV_ref, which contains the partial derivatives of slack power with respect to bus voltages. The function interprets changes in slack active power injections as indicative of changes in grid active power losses. KLU is used to factorize the sparse Jacobian matrix to solve for the loss factors.
Arguments
data::ACPowerFlowData: The data structure containing power flow information, including theloss_factorsmatrix.Jv::SparseMatrixCSC{Float64, Int32}: The sparse Jacobian matrix of the power flow system.time_step::Int: The time step index for which the loss factors are calculated.
PowerFlows._calculate_voltage_stability_factors — Method
calculate_voltage_stability_factors(data::ACPowerFlowData, J::ACPowerFlowJacobian, time_step::Integer)Calculate and store the voltage stability factors in the voltage_stability_factors matrix of the ACPowerFlowData structure for a given time step. The voltage stability factors are computed using the Jacobian matrix J in block format after a converged power flow calculation. The results are stored in the voltage_stability_factors matrix in the data instance. The factor for the grid as a whole (σ) is stored in the position of the REF bus. The values of the singular vector v indicate the sensitivity of the buses and are stored in the positions of the PQ buses. The values of v for PV buses are set to zero. The function uses the method described in "Fast calculation of a voltage stability index" by PA Lof et. al.
Arguments
data::ACPowerFlowData: The instance containing the grid model data.J::ACPowerFlowJacobian: The Jacobian matrix cache.time_step::Integer: The calculated time step.
PowerFlows._create_jacobian_matrix_structure — Method
_create_jacobian_matrix_structure(data::ACPowerFlowData, time_step::Int64) -> SparseMatrixCSC{Float64, Int32}Create the structure of the Jacobian matrix for an AC power flow problem.
Arguments
data::ACPowerFlowData: The power flow model.time_step::Int64: The specific time step for which the Jacobian matrix structure is created.
Returns
SparseMatrixCSC{Float64, Int32}: A sparse matrix with structural zeros representing the structure of the Jacobian matrix.
Description
This function initializes the structure of the Jacobian matrix for an AC power flow problem. The Jacobian matrix is used in power flow analysis to represent the partial derivatives of bus active and reactive power injections with respect to bus voltage magnitudes and angles.
Unlike some commonly used approaches where the Jacobian matrix is constructed as four submatrices, each grouping values for the four types of partial derivatives, this function groups the partial derivatives by bus. The structure is organized as groups of 4 values per bus.
This approach is more memory-efficient. Furthermore, this structure results in a more efficient factorization because the values are more likely to be grouped close to the diagonal. Refer to Electric Energy Systems: Analysis and Operation by Antonio Gomez-Exposito and Fernando L. Alvarado for more details.
The function initializes three arrays (rows, columns, and values) to store the row indices, column indices, and values of the non-zero elements of the Jacobian matrix, respectively.
For each bus in the system, the function iterates over its neighboring buses and determines the type of each neighboring bus (REF, PV, or PQ). Depending on the bus type, the function adds the appropriate entries to the Jacobian matrix structure.
- For
REFbuses, entries are added for local active and reactive power. - For
PVbuses, entries are added for active and reactive power with respect to angle, and for local reactive power. - For
PQbuses, entries are added for active and reactive power with respect to voltage magnitude and angle.
Example Structure
For a system with 3 buses where bus 1 is REF, bus 2 is PV, and bus 3 is PQ:
Let $\Delta P_j$, $\Delta Q_j$ be the active, reactive power balance at the $j$th bus. Let $P_j$ and $Q_j$ be the active and reactive power generated at the $j$th bus (REF and PV only). The state vector is $x = [P_1, Q_1, Q_2, \theta_2, V_3, \theta_3]$, and the residual vector is $F(x) = [\Delta P_1, \Delta Q_1, \Delta P_2, \Delta Q_2, \Delta P_3, \Delta Q_3]$.
The Jacobian matrix $J = \nabla F(x)$ has the structure:
\[J = \begin{bmatrix} \frac{\partial \vec{F}}{\partial P_1} & \frac{\partial \vec{F}}{\partial Q_1} & \frac{\partial \vec{F}}{\partial Q_2} & \frac{\partial \vec{F}}{\partial \theta_2} & \frac{\partial \vec{F}}{\partial V_3} & \frac{\partial \vec{F}}{\partial \theta_3} \end{bmatrix}\]
In reality, for large networks, this matrix would be sparse, and each 2×2 block would only be nonzero when there's a line between the respective buses.
Finally, the function constructs a sparse matrix from the collected indices and values and returns it.
PowerFlows._create_jacobian_matrix_structure_bus! — Method
Create the Jacobian matrix structure for a PV bus. Currently unused: we fill all four values even for PV buses with structiural zeros using the same function as for PQ buses.
PowerFlows._create_jacobian_matrix_structure_bus! — Method
Create the Jacobian matrix structure for a reference bus (REF). Currently unused: we fill all four values even for PV buses with structiural zeros using the same function as for PQ buses.
PowerFlows._create_jacobian_matrix_structure_bus! — Method
Create the Jacobian matrix structure for a PQ bus. Using this for all buses because a) for REF buses it doesn't matter if there are 2 values or 4 values - there are not many of them in the grid b) for PV buses we fill all four values because we can have a PV -> PQ transition and then we need to fill all four values
PowerFlows._create_jacobian_matrix_structure_lcc — Method
_create_jacobian_matrix_structure_lcc(
data::ACPowerFlowData,
rows::Vector{Int32},
columns::Vector{Int32},
values::Vector{Float64},
num_buses::Int
)Create the Jacobian matrix structure for LCC HVDC systems.
Description
The function iterates over each LCC system and adds the non-zero entries to the Jacobian matrix structure. The state vector for every LCC contains 4 variables: tap position and thyristor angle for both the rectifier and inverter sides. The indices of non-zero entries correspond to the positions of these variables in the extended state vector.
For an LCC system connecting bus $i$ (rectifier side) and bus $j$ (inverter side), the state variables are:
- $t_i$: tap position at rectifier
- $t_j$: tap position at inverter
- $\alpha_i$: thyristor angle at rectifier
- $\alpha_j$: thyristor angle at inverter
The residuals include:
- $F_{t_i}$: Active power balance at rectifier (controls $P_i$ to match setpoint)
- $F_{t_j}$: Total active power balance across LCC system
- $F_{\alpha_i}$: Rectifier thyristor angle constraint (maintains $\alpha_i$ at minimum)
- $F_{\alpha_j}$: Inverter thyristor angle constraint (maintains $\alpha_j$ at minimum)
Example Structure
For a system with 2 buses connected by one LCC where bus 1 is the rectifier side and bus 2 is the inverter side, the Jacobian matrix would have non-zero entries at positions like:
\[\begin{array}{c|cccccccc} & V_1 & \delta_1 & V_2 & \delta_2 & t_1 & t_2 & \alpha_1 & \alpha_2 \\ \hline P_1 & \frac{\partial P_1}{\partial V_1} & & & & \frac{\partial P_1}{\partial t_1} & & \frac{\partial P_1}{\partial \alpha_1} & \\ Q_1 & \frac{\partial Q_1}{\partial V_1} & & & & \frac{\partial Q_1}{\partial t_1} & & \frac{\partial Q_1}{\partial \alpha_1} & \\ P_2 & & & & & & & & \\ Q_2 & & & & & & & & \\ F_{t_1} & \frac{\partial F_{t_1}}{\partial V_1} & & & & \frac{\partial F_{t_1}}{\partial t_1} & & \frac{\partial F_{t_1}}{\partial \alpha_1} & \\ F_{t_2} & \frac{\partial F_{t_2}}{\partial V_1} & & \frac{\partial F_{t_2}}{\partial V_2} & & \frac{\partial F_{t_2}}{\partial t_1} & \frac{\partial F_{t_2}}{\partial t_2} & \frac{\partial F_{t_2}}{\partial \alpha_1} & \frac{\partial F_{t_2}}{\partial \alpha_2} \\ F_{\alpha_1} & & & & & & & \frac{\partial F_{\alpha_1}}{\partial \alpha_1} & \\ F_{\alpha_2} & & & & & & & & \frac{\partial F_{\alpha_2}}{\partial \alpha_2} \end{array}\]
This function sets up the indices of these non-zero entries in the sparse Jacobian matrix structure.
Arguments
data::ACPowerFlowData: The power flow data containing LCC system information.rows::Vector{Int32}: Vector to store row indices of non-zero Jacobian entries.columns::Vector{Int32}: Vector to store column indices of non-zero Jacobian entries.values::Vector{Float64}: Vector to store initial values of non-zero Jacobian entries.num_buses::Int: Total number of buses in the system.
PowerFlows._singular_value_decomposition — Method
_singular_value_decomposition(J::SparseMatrixCSC{Float64, Int32}, npvpq::Integer; tol::Float64 = 1e-9, max_iter::Integer = 100,)Estimate the smallest singular value σ and corresponding left and right singular vectors u and v of a sparse matrix G_s (a sub-matrix of J). This function uses an iterative method involving LU factorization of the Jacobian matrix to estimate the smallest singular value of G_s. The algorithm alternates between updating u and v, normalizing, and checking for convergence based on the change in the estimated singular value σ. The function uses the method described in Algorithm 3 of "Fast calculation of a voltage stability index" by PA Lof et. al.
Arguments
J::SparseMatrixCSC{Float64, Int32}: The sparse block-form Jacobian matrix.npvpq::Integer: Number of PV and PQ buses in J.
Keyword Arguments
tol::Float64=1e-9: Convergence tolerance for the iterative algorithm.max_iter::Integer=100: Maximum number of iterations.
Returns
σ::Float64: The estimated smallest singular value.left::Vector{Float64}: The estimated left singular vector (referred to asuin the cited paper).right::Vector{Float64}: The estimated right singular vector (referred to asvin the cited paper).
PowerFlows._update_jacobian_matrix_values! — Method
Used to update Jv based on the bus voltages, angles, etc. in data.