Internal

PowerFlows.ACPowerFlowJacobianType
struct ACPowerFlowJacobian

A 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 function Jf!.
source
PowerFlows.ACPowerFlowJacobianMethod
(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

J = ACPowerFlowJacobian(data, time_step)
J(time_step)  # Updates the Jacobian matrix Jv
source
PowerFlows.ACPowerFlowJacobianMethod
ACPowerFlowJacobian(data::ACPowerFlowData, time_step::Int64) -> ACPowerFlowJacobian

This 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.

Arguments

  • data::ACPowerFlowData: The data used for power flow calculations.
  • time_step::Int64: The time step for the calculations.

Returns

  • ACPowerFlowJacobian: An instance of ACPowerFlowJacobian.

#Example

J = ACPowerFlowJacobian(data, time_step)  # Creates an instance J of ACPowerFlowJacobian, with the Jacobian matrix stored internally as J.Jv initialized and updated.
J(time_step)  # Updates the Jacobian matrix stored internally in J (J.Jv) with the latest state of the `data` (`ACPowerFlowData` instance) and the provided time step.
J.Jv  # Access the Jacobian matrix stored internally in J.
source
PowerFlows.ACPowerFlowJacobianMethod
(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

J = ACPowerFlowJacobian(data, time_step)
Jv = SparseArrays.sparse(Float64[], Int32[], Int32[])
J(Jv, time_step)  # Updates the Jacobian matrix Jv and writes it to J
source
PowerFlows.ACPowerFlowResidualType
struct ACPowerFlowResidual

A 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.
source
PowerFlows.ACPowerFlowResidualMethod
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.

source
PowerFlows.ACPowerFlowResidualMethod
(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.
source
PowerFlows.ACPowerFlowResidualMethod
(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.
source
PowerFlows.KLULinSolveCacheType

A cached linear solver using KLU. Carefully written so as to minimize allocations: solve! and numeric_refactor! are completely non-allocating.

Fields:

  • K: the underlying KLU object.
  • reuse_symbolic::Bool: reuse the symbolic factorization. Defaults to true.
  • check_pattern::Bool: if true, numeric_refactor! verifies that the new

matrix has the same sparsity structure. Defaults to true. -rf_common, rf_symbolic, rf_numeric: internal usage. Stored to avoid allocations.

source
PowerFlows.LinearSolverCacheType

Abstract supertype for all cached linear solvers. Subtypes must implement: symbolic_factor!, symbolic_refactor!, numeric_refactor! (which doubles as numeric_factor!), and solve!.

source
PowerFlows.StateVectorCacheType

Cache for non-linear methods

Fields

-x::Vector{Float64}: the current state vector. -r::Vector{Float64}: the current residual. -Δx_nr::Vector{Float64}: the step under the Newton-Raphson method. The remainder of the fields are only used in the TrustRegionACPowerFlow: -r_predict::Vector{Float64}: the predicted residual at x+Δx_proposed, under a linear approximation: i.e J_x⋅(x+Δx_proposed). -Δx_proposed::Vector{Float64}: the suggested step Δx, selected among Δx_nr, Δx_cauchy, and the dogleg interpolation between the two. The first is chosen when x+Δx_nr is inside the trust region, the second when both x+Δx_cauchy and x+Δx_nr are outside the trust region, and the third when x+Δx_cauchy is inside and x+Δx_nr outside. The dogleg step selects the point where the line from x+Δx_cauchy to x+Δx_nr crosses the boundary of the trust region. -Δx_cauchy::Vector{Float64}: the step to the Cauchy point if the Cauchy point lies within the trust region, otherwise a step in that direction.

source
PowerFlows.A_plus_eq_BT_B!Method

Does A += B' * B, in a way that preserves the sparse structure of A, if possible. A workaround for the fact that Julia seems to run dropzeros!(A) automatically if I just do A .+= B' * B.

source
PowerFlows._calculate_loss_factorsMethod
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 injection 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 the loss_factors matrix.
  • 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.
source
PowerFlows._calculate_voltage_stability_factorsMethod
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 the following publication:

P.-A. Lof, T. Smed, G. Andersson, and D. J. Hill, "Fast calculation of a voltage stability index," in IEEE Transactions on Power Systems, vol. 7, no. 1, pp. 54-64, Feb. 1992, doi: 10.1109/59.141687.

Arguments

  • data::ACPowerFlowData: The instance containing the grid model data.
  • J::ACPowerFlowJacobian: The Jacobian matrix cache.
  • time_step::Integer: The calculated time step.
source
PowerFlows._create_jacobian_matrix_structureMethod
_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. Inputs are the grid model as an instance of ACPowerFlowData at a given time step.

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. See the example below for details.

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 REF buses, entries are added for local active and reactive power.
  • For PV buses, entries are added for active and reactive power with respect to angle, and for local reactive power.
  • For PQ buses, entries are added for active and reactive power with respect to voltage magnitude and angle.

For example, suppose we have a system with 3 buses: bus 1 is REF, bus 2 is PV, and bus 3 is PQ. Let ΔPⱼ, ΔQⱼ be the active, reactive power balance at the jth bus. Let Pⱼ and Qⱼ be the active and reactive power generated at the jth bus (REF and PV only). Then the state vector is [P₁, Q₁, Q₂, θ₂, V₃, θ₃], and the Jacobian matrix is

| ∂ΔP₁/∂P₁ | ∂ΔP₁/∂Q₁ | ∂ΔP₁/∂Q₂ | ∂ΔP₁/∂θ₂ | ∂ΔP₁/∂V₃ | ∂ΔP₁/∂θ₃ | | ∂ΔQ₁/∂P₁ | ∂ΔQ₁/∂Q₁ | ∂ΔQ₁/∂Q₂ | ∂ΔQ₁/∂θ₂ | ∂ΔQ₁/∂V₃ | ∂ΔQ₁/∂θ₃ | | ∂ΔP₂/∂P₁ | ∂ΔP₂/∂Q₁ | ∂ΔP₂/∂Q₂ | ∂ΔP₂/∂θ₂ | ∂ΔP₂/∂V₃ | ∂ΔP₂/∂θ₃ | | ∂ΔQ₂/∂P₁ | ∂ΔQ₂/∂Q₁ | ∂ΔQ₂/∂Q₂ | ∂ΔQ₂/∂θ₂ | ∂ΔQ₂/∂V₃ | ∂ΔQ₂/∂θ₃ | | ∂ΔP₃/∂P₁ | ∂ΔP₃/∂Q₁ | ∂ΔP₃/∂Q₂ | ∂ΔP₃/∂θ₂ | ∂ΔP₃/∂V₃ | ∂ΔP₃/∂θ₃ | | ∂ΔQ₃/∂P₁ | ∂ΔQ₃/∂Q₁ | ∂ΔQ₃/∂Q₂ | ∂ΔQ₃/∂θ₂ | ∂ΔQ₃/∂V₃ | ∂ΔQ₃/∂θ₃ |

In reality, for large networks, this matrix would be sparse, and each 2x2 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.

source
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

source
PowerFlows._dc_powerflow_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 powerflow on data::ACPowerFlowData for the given time_step, and writes the solution to data.bus_angles.

source
PowerFlows._dogleg!Method

Sets Δx_proposed equal to the Δx by which we should update x. Decides between the Cauchy step Δx_cauchy, Newton-Raphson step Δx_nr, and the dogleg interpolation between the two, based on which fall within the trust region.

source
PowerFlows._first_choice_gen_idMethod

Try to make an informative one or two character name for the load/generator/etc.

  • "generator-1234-AB" -> "AB"
  • "123CT7" -> "7"
  • "load1234" -> "34"
source
PowerFlows._fix_3w_transformer_ratingMethod

Setting a value of zero 0.0 when having a value greater than or equal to INFINITE_BOUND reverses the operation done in the PSY parsing side, according to PSSE Manual.

source
PowerFlows._make_gens_from_hvdcMethod

Create a synthetic generator (PSY.ThermalStandard) representing one end of a TwoTerminalGenericHVDCLine for export purposes. The generator is initialized with parameters reflecting the HVDC line's state.

Notes

- The generator's name is constructed as "<hvdc_line_name>_<suffix>".
- The `ext` field includes `"HVDC_END"` to indicate the end ("FR"/"TO").
source
PowerFlows._newton_powerflowMethod

Driver for the LevenbergMarquardtACPowerFlow method: sets up the data structures (e.g. residual), runs the powerflow method via calling _run_powerflow_method on them, then handles post-processing (e.g. loss factors).

source
PowerFlows._psse_bus_namesMethod

Given a vector of Sienna bus names, create a dictionary from Sienna bus name to PSS/E-compatible bus name. Guarantees determinism and minimal changes.

WRITTEN TO SPEC: PSS/E 33.3 POM 5.2.1 Bus Data

source
PowerFlows._psse_bus_numbersMethod

Given a vector of Sienna bus numbers, create a dictionary from Sienna bus number to PSS/E-compatible bus number. Assumes that the Sienna bus numbers are positive and unique. Guarantees determinism: if the input contains the same bus numbers in the same order, the output will. Guarantees minimal changes: that if an existing bus number is compliant, it will not be changed.

WRITTEN TO SPEC: PSS/E 33.3 POM 5.2.1 Bus Data

source
PowerFlows._psse_transformer_namesMethod

Given a vector of Sienna transformer names, create a dictionary from Sienna transformer name to PSS/E-compatible transformer name. Guarantees determinism and minimal changes.

WRITTEN TO SPEC: PSS/E 33.3/35.4 POM 5.2.1 Transformer Data

source
PowerFlows._run_powerflow_methodMethod

Runs the full NewtonRaphsonACPowerFlow.

Keyword arguments:

  • maxIterations::Int: maximum iterations. Default: 50.
  • tol::Float64: tolerance. The iterative search ends when norm(abs.(residual)) < tol. Default: 1.0e-9.
  • refinement_threshold::Float64: If the solution to J_x Δx = r satisfies norm(J_x Δx - r, 1)/norm(r, 1) > refinement_threshold, do iterative refinement to improve the accuracy. Default: 0.05.
  • refinement_eps::Float64: run iterative refinement on J_x Δx = r until norm(Δx_{i}-Δx_{i+1}, 1)/norm(r,1) < refinement_eps. Default: 1.0e-6
source
PowerFlows._run_powerflow_methodMethod

Runs the full TrustRegionNRMethod.

Keyword arguments:

  • maxIterations::Int: maximum iterations. Default: 50.
  • tol::Float64: tolerance. The iterative search ends when maximum(abs.(residual)) < tol. Default: 1.0e-9.
  • factor::Float64: the trust region starts out with radius factor*norm(x_0, 1), where x_0 is our initial guess, taken from data. Default: 1.0.
  • eta::Float64: improvement threshold. If the observed improvement in our residual exceeds eta times the predicted improvement, we accept the new x_i. Default: 0.0001.
source
PowerFlows._set_series_voltages_and_flows!Method

Calculate series voltages at buses removed in degree 2 reduction. Method: number the nodes in the series segment 0, 1, ..., n. Number the segments by their concluding node: 1, 2, ... n. The currents in the segments are given by: [y^iff y^ift; y^itf y^itt] * [V{i-1}; Vi] = [I{i-1, i}; I{i, i-1}] where I'm using upper indices to denote the segment number. There are no loads or generators at the internal nodes, so I{i, i+1} + I{i, i-1} = 0. Substitute the above expressions for the currents and group by Vi: y^i{tf} V{i-1} + (y{tt}^i + y{ff}^{i+1}) Vi + y{ft}^{i+1} V{i+1} = 0 For i = 1 and i = n-1, move the terms involving V0 and Vn [known] to the other side. This gives a tridiagonal system for x = [V1, ..., V{n-1}]: A * x = [-y^1{tf} * V0, 0, ..., 0, -y^{n}{ft} * Vn] where A has diagonal entries y{tt}^i + y{ff}^{i+1}, subdiagonal entries y{tf}^{i+1}, and superdiagonal entries y{ft}^i.

In the below, I use y11 instead of yff, y12 instead of yft, etc.

source
PowerFlows._simple_stepFunction

Does a single iteration of NewtonRaphsonACPowerFlow. Updates the r and x fields of the stateVector, and computes the Jacobian at the new x.

source
PowerFlows._singular_value_decompositionMethod
_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 in the following publication:

P.-A. Lof, T. Smed, G. Andersson, and D. J. Hill, "Fast calculation of a voltage stability index," in IEEE Transactions on Power Systems, vol. 7, no. 1, pp. 54-64, Feb. 1992, doi: 10.1109/59.141687.

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 as u in the cited paper).
  • right::Vector{Float64}: The estimated right singular vector (referred to as v in the cited paper).
source
PowerFlows._solve_Δx_nr!Method

Solve for the Newton-Raphson step, given the factorization object for J.Jv (if non-singular) or its stand-in (if singular).

source
PowerFlows._trust_region_stepMethod

Does a single iteration of the TrustRegionNRMethod: updates the x and r fields of the stateVector and computes the value of the Jacobian at the new x, if needed. Unlike _simple_step, this has a return value, the updated value of delta`.

source
PowerFlows._update_gens_from_hvdc!Method

Update the parameters of synthetic generators created from HVDC lines, so they reflect the current setpoints and limits of the HVDC devices in the system.

source
PowerFlows._update_hessian_matrix_values!Method

Sets Hv equal to F_1(x) H_{F_1}(x) + ...+ F_{2n}(x) H_{F_{2n}}(x), where Fk denotes the kth power balance equation and `H{F_k}its Hessian. This isn't the full Hessian of our function: it's only the terms in that come from the second derivatives of our power balance equations. (There's also aJ'*J` term.)

What's the sparse structure of that expression? It's split into 2x2 blocks, each corresponding to a pair of buses. The sparse structure of a block for a pair of buses connected by a branch is: | REF| PV | PQ –-+––+––+–– REF| | | | | | –-+––+––+–– PV | | | | | .| . . –-+––+––+–– PQ | | .| . . | | .| . . Diagonal blocks follow the same pattern as above (as if each bus is its own neighbor). Off-diagonal blocks for a pair of buses not connected by a branch are structurally zero.

source
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.
source
PowerFlows.block_J_indicesMethod
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:

| dPdθ | dPdV | | dQdθ | dQdV |

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.
source
PowerFlows.can_be_PVMethod

Return set of all bus numbers that can be PV: i.e. have an available generator, or certain voltage regulation devices.

source
PowerFlows.convert_emptyMethod

If val is empty, returns T(); if not, asserts that val isa T and returns val. Has nice type checker semantics.

Examples

convert_empty(Vector{String}, [])  # -> String[]
convert_empty(Vector{String}, ["a"])  # -> ["a"]
convert_empty(Vector{String}, [2])  # -> TypeError: in typeassert, expected Vector{String}, got a value of type Vector{Int64}
Base.return_types(Base.Fix1(convert_empty, Vector{String}))  # -> [Vector{String}]
source
PowerFlows.create_component_idsMethod

Given a vector of component names and a corresponding vector of container IDs (e.g., bus numbers), create unique-per-container PSS/E-compatible IDs, output a dictionary from (container ID, component name) to PSS/E-compatible component ID. The "singlesto1" flag detects components that are the only one on their bus and gives them the name "1".

source
PowerFlows.dc_powerflow_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.

source
PowerFlows.get_arc_namesMethod

Return the names of the arcs in the power flow data: those that correspond to branches in the system will get the branch names, others will get a placeholder name of the form from-to.

source
PowerFlows.get_branches_with_numbersMethod

Collects all AC branches (Line, MonitoredLine, DiscreteControlledACBranch) from the system, sorts them by their bus numbers, and returns a vector of tuples (branch, bus_numbers).

Arguments

  • exporter::PSSEExporter: The exporter containing the system.

Returns

  • Vector{Tuple{<:PSY.Branch, Tuple}}: Each tuple contains a branch and its associated bus numbers.
source
PowerFlows.make_power_flow_containerFunction

Create an appropriate PowerFlowContainer for the given PowerFlowEvaluationModel and initialize it from the given PSY.System.

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
  • time_steps::Int: number of time periods to consider (default is 1)
  • timestep_names::Vector{String}: names of the time periods defines by the argument "time_steps". Default value is String[].
source
PowerFlows.my_mul_mtMethod

Matrix multiplication A*x. Written this way because a VirtualPTDF matrix does not store all of its entries: instead, it calculates them (or retrieves them from cache), one element or one row at a time.

source
PowerFlows.numeric_refactor!Method

Frees numeric factorization stored by cache, if non-null. If cache.check_pattern is true and the sparse matrix structure of A doesn't match the cached one, throws an error. Finally, computes the numeric factorization of A and stores that to cache.

source
PowerFlows.partition_stateMethod

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.)

source
PowerFlows.set_values!Method
set_values!(mat::FixedStructureCHOLMOD, new_vals::AbstractVector{Float64})

In-place update of the numeric values in the CHOLMOD matrix.

source
PowerFlows.symbolic_factor!Method

Frees up the current symbolic and numeric factorizations stored by cache, if non-null. Then computes the symbolic factorization of A and stores that to cache.

source
PowerFlows.symbolic_refactor!Method

Symbolic refactor. Behavior depends on the values of cache.reuse_symbol and cache.check_pattern. There are 3 cases:

  • !reuse_symbol: always refactor. Just calls symbolic_factor(cache, A).
  • reuse_symbol && check_pattern: checks if the symbolic structure of A matches the cached one, and throws an error if it doesn't. This is to prevent bad input: we expected the structure to be the same, but it isn't.
  • reuse_symbol && !check pattern: do nothing. Assume the structure of A matches the cached one.
source
PowerFlows.update_system!Method

Modify the values in the given System to correspond to the given PowerFlowData such that if a new PowerFlowData is constructed from the resulting system it is the same as data. See also write_powerflow_solution!. NOTE that this assumes that data was initialized from sys and then solved with no further modifications.

source
PowerFlows.write_to_buffers!Method

WRITTEN TO SPEC: PSS/E 33.3/35.4 POM 5.2.1 Bus Data. Sienna voltage limits treated as PSS/E normal voltage limits; PSSE emergency voltage limits left as default.

source