API Reference

Sensitivity Interface

PowerDiff.calc_sensitivityFunction
calc_sensitivity(state, operand::Symbol, parameter::Symbol) → Sensitivity{T}

Compute sensitivity of operand with respect to parameter.

Returns a Sensitivity{T} result that:

  • Acts like a matrix (implements AbstractMatrix interface)
  • Has symbol fields for formulation, operand, and parameter
  • Includes bidirectional index mappings for element IDs

Invalid combinations throw ArgumentError.

Operand Symbols

  • :va: Voltage phase angles (DC PF, DC OPF, AC PF, AC OPF)
  • :f: Branch active power flows (DC PF, DC OPF, AC PF)
  • :pg / :g: Generator active power (DC OPF, AC OPF)
  • :psh: Load shedding (DC OPF only)
  • :qg: Generator reactive power (AC OPF)
  • :lmp: Locational marginal prices (DC OPF, AC OPF)
  • :qlmp: Reactive power locational marginal prices (AC OPF)
  • :vm: Voltage magnitude (AC PF, AC OPF)
  • :im: Current magnitude (AC PF)
  • :v: Complex voltage phasor (AC PF) — returns ComplexF64 elements
  • :p: Active power injection (AC PF, as operand for Jacobian blocks)
  • :q: Reactive power injection (AC PF, as operand for Jacobian blocks)

Parameter Symbols

  • :d / :pd: Active demand (DC PF, DC OPF, AC OPF; AC PF via transform)
  • :qd: Reactive demand (AC OPF; AC PF via transform)
  • :sw: Switching states
  • :cq, :cl: Cost coefficients (DC OPF, AC OPF)
  • :fmax: Flow limits (DC OPF, AC OPF)
  • :b: Branch susceptances (DC PF, DC OPF, AC PF)
  • :g: Branch conductances (AC PF)
  • :p: Active power injection (AC PF)
  • :q: Reactive power injection (AC PF)
  • :va: Voltage phase angle (AC PF, as parameter for Jacobian blocks)
  • :vm: Voltage magnitude (AC PF, as parameter for Jacobian blocks)

Examples

# DC Power Flow
pf_state = DCPowerFlowState(net, demand)
sens = calc_sensitivity(pf_state, :va, :d)
sens.formulation  # :dcpf
sens.operand      # :va
sens.parameter    # :d

# DC OPF
prob = DCOPFProblem(net, demand)
solve!(prob)
sens = calc_sensitivity(prob, :lmp, :d)

# AC Power Flow — voltage/current sensitivities
sens = calc_sensitivity(ac_state, :vm, :p)
topo = calc_sensitivity(ac_state, :vm, :g)

# AC Power Flow — Jacobian blocks
J1 = calc_sensitivity(ac_state, :p, :va)   # ∂P/∂θ
J2 = calc_sensitivity(ac_state, :p, :vm)   # ∂P/∂|V|

# AC Power Flow — demand transform (∂/∂d = -∂/∂p since p_net = pg - pd)
dvm_dd = calc_sensitivity(ac_state, :vm, :d)

# AC OPF — all parameter types supported
dlmp_dd = calc_sensitivity(ac_prob, :lmp, :d)      # dLMP/dd
dpg_dcq = calc_sensitivity(ac_prob, :pg, :cq)      # dpg/dcq
dvm_dfmax = calc_sensitivity(ac_prob, :vm, :fmax)   # d|V|/dfmax
source
PowerDiff.calc_sensitivity_columnFunction
calc_sensitivity_column(state, operand::Symbol, parameter::Symbol, col_id::Int) → Vector

Compute a single column of the sensitivity matrix — the sensitivity of operand with respect to a single element of parameter identified by col_id.

col_id is an element ID (bus/branch/gen ID), consistent with Sensitivity.col_to_id. Returns a dense vector of length matching the operand dimension.

For OPF problems, this avoids materializing the full N×N sensitivity matrix by solving a single KKT backsubstitution (O(nnz) instead of O(nnz × n)).

Examples

# Sensitivity of all LMPs to demand at bus 3
col = calc_sensitivity_column(prob, :lmp, :d, 3)

# Compare with full matrix
S = calc_sensitivity(prob, :lmp, :d)
col ≈ Matrix(S)[:, S.id_to_col[3]]  # true
source
PowerDiff.SensitivityType
Sensitivity{T} <: AbstractMatrix{T}

Sensitivity matrix with symbol metadata for formulation, operand, and parameter.

Wraps a matrix with bidirectional index mappings between internal indices and original element IDs. Implements the full AbstractMatrix interface.

Type Parameter

  • T <: Number: Element type (Float64, ComplexF64)

Fields

  • matrix::Matrix{T}: The sensitivity matrix data
  • formulation::Symbol: Formulation tag (:dcpf, :dcopf, :acpf, :acopf)
  • operand::Symbol: Operand tag (:va, :vm, :pg, :qg, :f, :lmp, :psh, :im, :v)
  • parameter::Symbol: Parameter tag (:d, :sw, :cq, :cl, :fmax, :b, :g, :p, :q, :va, :vm, :qd)
  • row_to_id::Vector{Int}: Maps internal row index → external element ID
  • id_to_row::Dict{Int,Int}: Maps external element ID → internal row index
  • col_to_id::Vector{Int}: Maps internal col index → external element ID
  • id_to_col::Dict{Int,Int}: Maps external element ID → internal col index

Examples

sens = calc_sensitivity(prob, :lmp, :d)
sens.formulation  # :dcopf
sens.operand      # :lmp
sens.parameter    # :d
size(sens)        # (n, n)
sens[2, 3]        # Access element
Matrix(sens)      # Get raw matrix
source

ID Mapping

PowerDiff.IDMappingType
IDMapping

Bidirectional mapping between original PowerModels element IDs and sequential 1-based indices used for internal computation.

Constructed from a PM.build_ref dictionary or as an identity mapping for direct/programmatic constructors.

Fields

  • bus_ids::Vector{Int}: Sorted original bus IDs; index position = sequential index
  • branch_ids::Vector{Int}: Sorted original branch IDs
  • gen_ids::Vector{Int}: Sorted original generator IDs
  • load_ids::Vector{Int}: Sorted original load IDs
  • shunt_ids::Vector{Int}: Sorted original shunt IDs
  • bus_to_idx::Dict{Int,Int}: Original bus ID → sequential index
  • branch_to_idx::Dict{Int,Int}: Original branch ID → sequential index
  • gen_to_idx::Dict{Int,Int}: Original generator ID → sequential index
  • load_to_idx::Dict{Int,Int}: Original load ID → sequential index
  • shunt_to_idx::Dict{Int,Int}: Original shunt ID → sequential index
source

JVP / VJP

PowerDiff.jvpFunction
jvp(S::Sensitivity, δp::AbstractDict{Int,<:Number}) → Dict{Int,T}

Jacobian-vector product with ID-aware input and output. Computes S * δp where δp is keyed by original parameter element IDs. Returns a Dict keyed by original operand element IDs.

source
jvp(S::Sensitivity, δp::AbstractVector) → Dict{Int,T}

Jacobian-vector product with sequential input and ID-aware output.

source
jvp(state, operand::Symbol, parameter::Symbol, tang::AbstractVector) → Vector

Efficient Jacobian-vector product (∂operand/∂parameter) · tang without materializing the full sensitivity matrix.

Uses a single KKT forward-solve (O(nnz)) instead of building the full O(n²) sensitivity matrix. For performance-critical loops, use jvp! with pre-allocated buffers.

Examples

prob = DCOPFProblem(net, d); solve!(prob)
δd = randn(net.n)
δlmp = jvp(prob, :lmp, :d, δd)   # ∂LMP/∂d · δd, length n
source
jvp(state, operand::Symbol, parameter::Symbol, tang::AbstractDict{Int}) → Dict{Int,Float64}

ID-aware JVP. Input keyed by parameter element IDs, output keyed by operand element IDs.

Note: allocates a Dict for the output. For performance-critical loops, prefer the AbstractVector interface.

source
PowerDiff.jvp!Function
jvp!(out, prob::DCOPFProblem, operand::Symbol, parameter::Symbol,
     tang::AbstractVector; work=nothing) → out

In-place JVP for performance-critical loops (e.g., bilevel optimization).

out must be pre-allocated to the operand dimension. Pass work (a Vector{Float64} of length kkt_dims(prob)) to override the internal workspace (e.g., for thread-local buffers); if omitted, a workspace is lazily allocated in prob.cache and reused across calls.

Examples

prob = DCOPFProblem(net, d); solve!(prob)
out = zeros(net.n)

# Simple (workspace managed by cache):
jvp!(out, prob, :lmp, :d, tang)

# Thread-local override:
work = zeros(kkt_dims(prob))
jvp!(out, prob, :lmp, :d, tang; work=work)
source
PowerDiff.vjpFunction
vjp(S::Sensitivity, δy::AbstractDict{Int,<:Number}) → Dict{Int,T}

Vector-Jacobian product with ID-aware input and output. Computes S' * δy where δy is keyed by original operand element IDs. Returns a Dict keyed by original parameter element IDs.

source
vjp(S::Sensitivity, δy::AbstractVector) → Dict{Int,T}

Vector-Jacobian product with sequential input and ID-aware output.

source
vjp(state, operand::Symbol, parameter::Symbol, adj::AbstractVector) → Vector

Efficient vector-Jacobian product (∂operand/∂parameter)ᵀ · adj without materializing the full sensitivity matrix.

Uses a single KKT transpose-solve (O(nnz)) instead of building the full O(n²) sensitivity matrix. For performance-critical loops, use vjp! with pre-allocated buffers.

Examples

prob = DCOPFProblem(net, d); solve!(prob)
w = randn(net.n)
grad_d = vjp(prob, :lmp, :d, w)   # ∂LMP/∂dᵀ · w, length n
source
vjp(state, operand::Symbol, parameter::Symbol, adj::AbstractDict{Int}) → Dict{Int,Float64}

ID-aware VJP. Input keyed by operand element IDs, output keyed by parameter element IDs.

Note: allocates a Dict for the output. For performance-critical loops, prefer the AbstractVector interface.

source
PowerDiff.vjp!Function
vjp!(out, prob::DCOPFProblem, operand::Symbol, parameter::Symbol,
     adj::AbstractVector; work=nothing) → out

In-place VJP for performance-critical loops (e.g., bilevel optimization).

out must be pre-allocated to the parameter dimension. Pass work (a Vector{Float64} of length kkt_dims(prob)) to override the internal workspace (e.g., for thread-local buffers); if omitted, a workspace is lazily allocated in prob.cache and reused across calls.

Examples

prob = DCOPFProblem(net, d); solve!(prob)
out = zeros(net.n)

# Simple (workspace managed by cache):
vjp!(out, prob, :lmp, :d, adj)

# Thread-local override:
work = zeros(kkt_dims(prob))
vjp!(out, prob, :lmp, :d, adj; work=work)
source
PowerDiff.kkt_dimsFunction
kkt_dims(prob::AbstractOPFProblem) → Int

Total dimension of the KKT system for the given OPF problem.

source
PowerDiff.dict_to_vecFunction
dict_to_vec(S::Sensitivity, d::AbstractDict{Int,<:Number}, axis::Symbol) → Vector

Convert a Dict keyed by original element IDs to a dense vector aligned with the sensitivity matrix. Missing keys are treated as zero. Unknown IDs throw ArgumentError.

axis is :col (parameter space) or :row (operand space).

source
PowerDiff.vec_to_dictFunction
vec_to_dict(S::Sensitivity{T}, v::AbstractVector, axis::Symbol) → Dict{Int,T}

Convert a dense vector to a Dict keyed by original element IDs.

axis is :col (parameter space) or :row (operand space).

source

Introspection

PowerDiff.operand_symbolsFunction
operand_symbols(state) → Vector{Symbol}

Return the valid operand symbols for calc_sensitivity on the given state or problem. Includes symbols available via parameter transforms.

Examples

operand_symbols(pf_state)  # [:va, :f]
operand_symbols(prob)       # [:va, :pg, :f, :psh, :lmp]  (DCOPFProblem)
source
PowerDiff.parameter_symbolsFunction
parameter_symbols(state) → Vector{Symbol}

Return the valid parameter symbols for calc_sensitivity on the given state or problem. Includes symbols available via parameter transforms.

Examples

parameter_symbols(pf_state)  # [:d, :sw, :b]
parameter_symbols(prob)       # [:d, :sw, :cq, :cl, :fmax, :b]  (DCOPFProblem)
source

DC Types

PowerDiff.DCNetworkType
DCNetwork <: AbstractPowerNetwork

DC network data for B-theta OPF formulation. Uses susceptance-weighted Laplacian B = A' * Diagonal(-b .* sw) * A which preserves graphical structure for topology sensitivity analysis.

Fields

  • n, m, k: Number of buses, branches, and generators
  • A: Branch-bus incidence matrix (m x n)
  • G_inc: Generator-bus incidence matrix (n x k)
  • b: Branch susceptances (imaginary part of 1/z)
  • sw: Branch switching states (1 = closed, 0 = open)
  • fmax, gmax, gmin: Flow and generation limits
  • angmax, angmin: Phase angle difference limits
  • cq, cl: Quadratic and linear generation cost coefficients
  • c_shed: Load-shedding cost per bus (penalty for involuntary load curtailment)
  • ref_bus: Reference bus index (phase angle = 0)
  • tau: Regularization parameter for strong convexity
  • id_map: Bidirectional mapping between original and sequential element IDs
  • ref: PowerModels build_ref dictionary (nothing for programmatic constructors)
source
PowerDiff.DCPowerFlowStateType
DCPowerFlowState{F} <: AbstractPowerFlowState

DC power flow solution (phase angles from reduced-Laplacian solve, no optimization). Supports both generation and demand for flexible sensitivity analysis.

Unlike DCOPFSolution, this represents a simple power flow solution θ_r = B_r \ p_r where B_r is the susceptance matrix with the reference bus row and column deleted (invertible for a connected network), without optimal dispatch or constraint handling.

Fields

  • net: DCNetwork data
  • va: Phase angles (rad), with va[ref_bus] = 0
  • p: Net injection vector (p = pg - d)
  • pg: Generation vector
  • d: Demand vector
  • f: Branch flows (computed from va)
  • B_r_factor: Factorization of B[non_ref, non_ref] (Cholesky for inductive networks, LU fallback)
  • non_ref: Indices of non-reference buses
source
PowerDiff.DCOPFProblemType
DCOPFProblem <: AbstractOPFProblem

B-θ formulation of DC OPF wrapped around a JuMP model.

Fields

  • model: JuMP Model
  • network: DCNetwork data
  • va, pg, f, psh: Variable references for phase angles, generation, flows, load shedding
  • d: Demand parameter (can be updated for sensitivity analysis)
  • cons: Named tuple of constraint references
  • cache: Mutable sensitivity cache for avoiding redundant KKT solves
  • _optimizer: Optimizer factory for model rebuilds (internal)
  • _silent: Whether to suppress solver output (internal)
source
PowerDiff.DCOPFSolutionType
DCOPFSolution <: AbstractOPFSolution

Solution container for DC OPF problem, storing both primal and dual variables.

Fields

  • va: Phase angles at each bus
  • pg: Generator outputs
  • f: Line flows
  • psh: Load shedding at each bus
  • nu_bal: Power balance dual variables (nodal, used for LMP computation)
  • nu_flow: Flow definition dual variables
  • lam_ub, lam_lb: Line flow upper/lower bound duals
  • rho_ub, rho_lb: Generator upper/lower bound duals
  • mu_lb, mu_ub: Load shedding lower/upper bound duals
  • gamma_lb, gamma_ub: Phase angle difference lower/upper bound duals
  • objective: Optimal objective value
  • B_r_factor: Cached factorization of reduced susceptance matrix B[non_ref, non_ref]
source
PowerDiff.DCSensitivityCacheType
DCSensitivityCache

Mutable cache for DC OPF sensitivity data to avoid redundant KKT solves.

DC OPF supports 6 parameter types (:d, :sw, :cl, :cq, :fmax, :b), each producing a separate dz_d* full-derivative matrix. All share one KKT LU factorization (kkt_factor), so after the first parameter type is queried the factorization is reused for subsequent parameters. Different operand queries (e.g. :va vs :pg vs :lmp) for the same parameter type all extract rows from the same cached dz_d* matrix — no recomputation needed.

ACSensitivityCache follows the same design: one KKT factorization shared across 6 parameter types (:sw, :d, :qd, :cq, :cl, :fmax), each producing a separate cached dz_d* matrix. Power flow states (DCPowerFlowState, ACPowerFlowState) have no cache at all because their sensitivities are cheap direct algebra (reduced-Laplacian factorization or Jacobian factorization is precomputed at construction time).

Fields

  • solution: Cached DCOPFSolution (or nothing if not yet solved)
  • kkt_factor: Cached LU factorization of KKT Jacobian (or nothing)
  • dz_dd: Full KKT derivative w.r.t. demand (or nothing)
  • dz_dcl: Full KKT derivative w.r.t. linear cost (or nothing)
  • dz_dcq: Full KKT derivative w.r.t. quadratic cost (or nothing)
  • dz_dsw: Full KKT derivative w.r.t. switching (or nothing)
  • dz_dfmax: Full KKT derivative w.r.t. flow limits (or nothing)
  • dz_db: Full KKT derivative w.r.t. susceptances (or nothing)
  • b_r_factor: Cached reduced susceptance factorization (topology-dependent, survives demand changes)
  • work: Scratch workspace for VJP/JVP KKT solves (lazily allocated on first call, survives invalidation since its size depends only on (n, m, k) which are fixed at construction)
source

DC Functions

PowerDiff.solve!Function
solve!(prob::DCOPFProblem)

Solve the DC OPF problem and return a DCOPFSolution.

Invalidates the sensitivity cache since the solution may have changed.

Returns

DCOPFSolution containing optimal primal and dual variables.

Throws

Error if optimization does not converge to optimal/locally optimal solution.

source
solve!(prob::ACOPFProblem)

Solve the AC OPF problem and return an ACOPFSolution.

Invalidates the sensitivity cache since the solution may have changed.

Returns

ACOPFSolution containing optimal primal and dual variables.

Throws

Error if optimization does not converge to optimal/locally optimal solution.

source
PowerDiff.update_demand!Function
update_demand!(prob::DCOPFProblem, d::AbstractVector)

Update the demand parameter in the DC OPF problem.

This modifies the RHS of power balance and shedding upper-bound constraints for re-solving with new demand. Invalidates the sensitivity cache since parameters have changed.

source
PowerDiff.update_switching!Function
update_switching!(prob::DCOPFProblem, s::AbstractVector)

Update the network switching state, invalidate the sensitivity cache, and rebuild the JuMP model so that solve!(prob) uses the new switching state.

Arguments

  • prob: DCOPFProblem to update
  • s: New switching state vector (length m), values in [0,1]
source
update_switching!(prob::ACOPFProblem, sw::AbstractVector)

Update the network switching state, invalidate the sensitivity cache, and rebuild the JuMP model so that solve!(prob) uses the new switching state.

Arguments

  • prob: ACOPFProblem to update
  • sw: New switching state vector (length m), values in [0,1]
source
PowerDiff.calc_demand_vectorFunction
calc_demand_vector(net::Dict)

Extract demand vector from PowerModels network dictionary.

Works with both basic and non-basic networks. For non-basic networks, uses PM.build_ref to resolve load-bus mappings and returns a vector in sequential bus order (matching IDMapping from DCNetwork(net)).

source
calc_demand_vector(network::DCNetwork)

Extract demand vector from a DCNetwork that was constructed from a PowerModels dict.

Uses the stored ref to avoid redundant build_ref calls.

source
PowerDiff.calc_susceptance_matrixFunction
calc_susceptance_matrix(network::DCNetwork)

Compute the susceptance-weighted Laplacian: B = A' * Diagonal(-b .* sw) * A.

Sign convention: b stores Im(1/z) which is negative for inductive branches. The negation -b produces positive edge weights, making B positive semidefinite. This is the negative of PowerModels' calc_susceptance_matrix (which uses the standard bus susceptance matrix convention with negative diagonal).

DC power flow: B * θ = p (net injection). Branch flows: f = Diag(-b .* sw) * A * θ.

source
PowerDiff.calc_lmpFunction
calc_lmp(sol::DCOPFSolution, net::DCNetwork)

Compute Locational Marginal Prices from DC OPF solution.

The LMP at bus i is the marginal cost of serving an additional unit of demand at that bus. In the B-θ formulation, this equals the power balance dual ν_bal[i].

Returns

Vector of LMPs (length n), one per bus.

Example

sol = solve!(prob)
lmps = calc_lmp(sol, prob.network)
source
calc_lmp(prob::DCOPFProblem)

Solve the problem (if needed) and compute LMPs.

source
calc_lmp(sol::ACOPFSolution, prob::ACOPFProblem)

Compute Locational Marginal Prices from AC OPF solution.

The LMP at bus i is the marginal cost of serving an additional unit of active demand at that bus: LMPi = ∂f*/∂Pdi = -νpbali.

The sign negation arises because JuMP's dual ν_p_bal is negative at optimum for the standard power balance formulation (see sign derivation in file header).

Returns

Vector of LMPs (length n), one per bus.

source
calc_lmp(prob::ACOPFProblem)

Solve the AC OPF problem (if needed) and compute LMPs.

source
PowerDiff.calc_qlmpFunction
calc_qlmp(sol::ACOPFSolution, prob::ACOPFProblem)

Compute reactive power Locational Marginal Prices from AC OPF solution.

The QLMP at bus i is the marginal cost of serving an additional unit of reactive demand at that bus: QLMPi = ∂f*/∂Qdi = -νqbali.

Returns

Vector of QLMPs (length n), one per bus.

source
calc_qlmp(prob::ACOPFProblem)

Solve the AC OPF problem (if needed) and compute reactive power LMPs.

source
PowerDiff.calc_congestion_componentFunction
calc_congestion_component(sol::DCOPFSolution, net::DCNetwork)

Extract the congestion component of LMPs for analysis.

From the θ-stationarity KKT condition: B' * νbal + (WA)' * νflow + eref * ηref + A'*(γub - γlb) = 0

The congestion RHS includes both flow limit duals and angle difference duals: congestion[nonref] = Br \ (A' W (λub - λlb) + A'(γub - γlb))[non_ref]

The congestion component captures price differentiation due to binding flow and angle constraints, with the reference bus congestion component equal to zero.

Returns

Vector (length n) of congestion contributions to each bus's LMP.

source
PowerDiff.calc_energy_componentFunction
calc_energy_component(sol::DCOPFSolution, net::DCNetwork)

Extract the energy (non-congestion) component of LMPs for analysis.

This is the uniform price component: energy = ν_bal - congestion For a connected network, this should be approximately constant across all buses.

Returns

Vector (length n) of energy contributions to each bus's LMP.

source
PowerDiff.calc_generation_participation_factorsFunction
calc_generation_participation_factors(prob::DCOPFProblem) → Sensitivity{Float64}

Compute generation participation factors from demand sensitivity.

The participation factor for generator i at bus j is dgi/ddj, representing how much generator i output changes when demand at bus j increases by 1 MW.

Returns

Sensitivity{Float64} (k × n) with formulation=:dcopf, operand=:pg, parameter=:d.

source
PowerDiff.calc_ptdf_from_sensitivityFunction
calc_ptdf_from_sensitivity(prob::DCOPFProblem) → Sensitivity{Float64}

Compute Power Transfer Distribution Factors from flow sensitivity.

PTDF[e, j] = dfe/ddj represents how much flow on line e changes when power is injected at bus j (and withdrawn at the slack bus).

For the B-theta formulation, this is directly available from the sensitivity analysis.

Returns

Sensitivity{Float64} (m × n) with formulation=:dcopf, operand=:f, parameter=:d.

source
PowerDiff.ptdf_matrixFunction
ptdf_matrix(state::DCPowerFlowState) → Matrix{Float64}

Return the standard PTDF matrix (∂f/∂p) from a DC power flow state.

PD's calc_sensitivity(state, :f, :d) computes ∂f/∂d = -PTDF because p = g - d ⟹ ∂p/∂d = -I. This function negates to recover the standard PTDF sign convention: PTDF = ∂f/∂p.

source
PowerDiff.invalidate!Function
invalidate!(cache::DCSensitivityCache)

Clear cached sensitivity data that depends on the current solution. The b_r_factor field is preserved because it depends only on network topology (susceptances and switching), not on demand or the optimization solution.

source
invalidate!(cache::ACSensitivityCache)

Clear all cached AC sensitivity data. Called when problem parameters change.

source

AC Types

PowerDiff.ACNetworkType
ACNetwork <: AbstractPowerNetwork

AC network data with vectorized admittance representation.

Provides a unified interface for AC power flow and sensitivity analysis, analogous to DCNetwork for DC formulations. Uses edge-based conductance and susceptance vectors for differentiable admittance matrix construction.

The admittance matrix is reconstructed as: Y = A' * Diag(g + jb) * A + Diag(g_shunt + jb_shunt)

For switching-aware formulation: Y(sw) = A' * Diag((g + jb) . sw) * A + Diag(gshunt + j*bshunt)

Fields

  • n: Number of buses
  • m: Number of branches
  • A: Branch-bus incidence matrix (m × n)
  • incidences: Edge list [(i,j), ...] for each branch (sequential indices)
  • g: Branch conductances
  • b: Branch susceptances (note: typically negative for inductive lines)
  • g_shunt: Shunt conductances per bus (from shunts + line charging)
  • b_shunt: Shunt susceptances per bus
  • sw: Branch switching states ∈ [0,1]^m
  • is_switchable: Which branches can be switched
  • idx_slack: Slack bus index (sequential)
  • vm_min, vm_max: Voltage magnitude limits per bus
  • i_max: Branch current magnitude limits
  • id_map: Bidirectional mapping between original and sequential element IDs
  • ref: PowerModels build_ref dictionary (nothing for Y-matrix constructors)
source
PowerDiff.ACPowerFlowStateType
ACPowerFlowState <: AbstractPowerFlowState

AC power flow solution with full injection tracking.

Provides a common interface for AC sensitivity computations, analogous to DCPowerFlowState for DC power flow. Can be constructed from a PowerModels network or from raw voltage/admittance data.

Fields

  • net: ACNetwork reference (optional, provides access to edge-level data)
  • v: Complex voltage phasors at all buses
  • Y: Bus admittance matrix
  • p: Net real power injection (p = pg - pd)
  • q: Net reactive power injection (q = qg - qd)
  • pg: Real power generation per bus
  • pd: Real power demand per bus
  • qg: Reactive power generation per bus
  • qd: Reactive power demand per bus
  • branch_data: Branch dictionary with sequential indices (required for :im sensitivity)
  • idx_slack: Index of the slack (reference) bus
  • n: Number of buses
  • m: Number of branches

Constructors

  • ACPowerFlowState(v, Y; ...): From voltage phasors and admittance matrix
  • ACPowerFlowState(net::ACNetwork, v; ...): From ACNetwork and voltage solution
  • ACPowerFlowState(pm_net::Dict): From solved PowerModels network
source
PowerDiff.ACOPFProblemType
ACOPFProblem <: AbstractOPFProblem

Polar coordinate AC OPF wrapped around a JuMP model.

Fields

  • model: JuMP Model
  • network: ACNetwork data
  • va, vm: Variable references for voltage angles and magnitudes
  • pg, qg: Variable references for active and reactive generation
  • p, q: Dict of branch flow variable references (keyed by arc tuple)
  • cons: Named tuple of constraint references
  • ref: PowerModels-style reference dictionary (sequential keys)
  • gen_buses: Generator bus indices (maps generator index to bus index)
  • n_gen: Number of generators
  • cache: ACSensitivityCache for caching KKT derivatives
  • _optimizer: Optimizer factory for model rebuilds (internal)
  • _silent: Whether to suppress solver output (internal)
source
PowerDiff.ACOPFSolutionType
ACOPFSolution <: AbstractOPFSolution

Solution container for AC OPF problem, storing primal and dual variables.

Fields

Primal Variables

  • va: Voltage angles at each bus (n)
  • vm: Voltage magnitudes at each bus (n)
  • pg: Active power generation (k)
  • qg: Reactive power generation (k)
  • p: Active branch flows Dict{(l,i,j) => Float64}
  • q: Reactive branch flows Dict{(l,i,j) => Float64}

Dual Variables (Equality Constraints)

  • nu_p_bal: Active power balance duals (n) - used for LMP
  • nu_q_bal: Reactive power balance duals (n)
  • nu_ref_bus: Reference bus constraint duals (n_ref, usually 1)
  • nu_p_fr, nu_p_to: Active flow definition duals (m each)
  • nu_q_fr, nu_q_to: Reactive flow definition duals (m each)

Dual Variables (Inequality Constraints)

  • lam_thermal_fr, lam_thermal_to: Thermal limit duals (m each)
  • lam_angle_lb, lam_angle_ub: Angle difference limit duals (m each)
  • mu_vm_lb, mu_vm_ub: Voltage magnitude bound duals (n each)
  • rho_pg_lb, rho_pg_ub: Active gen bound duals (k each)
  • rho_qg_lb, rho_qg_ub: Reactive gen bound duals (k each)
  • sig_p_fr_lb, sig_p_fr_ub: From-side active flow bound duals (m each)
  • sig_q_fr_lb, sig_q_fr_ub: From-side reactive flow bound duals (m each)
  • sig_p_to_lb, sig_p_to_ub: To-side active flow bound duals (m each)
  • sig_q_to_lb, sig_q_to_ub: To-side reactive flow bound duals (m each)

Objective

  • objective: Optimal objective value
source
PowerDiff.ACSensitivityCacheType
ACSensitivityCache

Mutable cache for storing computed AC OPF sensitivity data to avoid redundant KKT solves and ForwardDiff Jacobian evaluations.

AC OPF supports 6 parameter types (:sw, :d, :qd, :cq, :cl, :fmax), each producing a separate dz_d* full-derivative matrix. All share one KKT LU factorization (kkt_factor), so after the first parameter type is queried the factorization is reused for subsequent parameters. Different operand queries (e.g. :va vs :pg vs :lmp) for the same parameter type all extract rows from the same cached dz_d* matrix — no recomputation needed.

Fields

  • solution: Cached ACOPFSolution (or nothing if not yet solved)
  • kkt_factor: Cached LU factorization of KKT Jacobian (or nothing)
  • kkt_constants: Cached KKT constants NamedTuple (or nothing)
  • dz_dsw: Full KKT derivative w.r.t. switching (or nothing)
  • dz_dd: Full KKT derivative w.r.t. active demand (or nothing)
  • dz_dqd: Full KKT derivative w.r.t. reactive demand (or nothing)
  • dz_dcq: Full KKT derivative w.r.t. quadratic cost (or nothing)
  • dz_dcl: Full KKT derivative w.r.t. linear cost (or nothing)
  • dz_dfmax: Full KKT derivative w.r.t. flow limits (or nothing)
source

AC Functions

PowerDiff.admittance_matrixFunction
admittance_matrix(net::ACNetwork) → SparseMatrixCSC{ComplexF64}

Reconstruct the bus admittance matrix Y from vectorized representation.

Y = A' * Diag(g + j*b) * A + Diag(g_shunt + j*b_shunt)
source
admittance_matrix(net::ACNetwork, sw::AbstractVector) → SparseMatrixCSC{ComplexF64}

Reconstruct admittance matrix with switching states.

Y(sw) = A' * Diag((g + j*b) .* sw) * A + Diag(g_shunt + j*b_shunt)
source
PowerDiff.branch_currentFunction
branch_current(net::ACNetwork, v::AbstractVector{<:Complex}) → Vector{ComplexF64}

Complex branch currents: I_branch = Diag(y) * A * v

source
PowerDiff.branch_powerFunction
branch_power(net::ACNetwork, v::AbstractVector{<:Complex}) → Vector{ComplexF64}

Complex branch power flows: Sbranch = diag(A*v) * conj(Ibranch)

source
PowerDiff.calc_power_flow_jacobianFunction
calc_power_flow_jacobian(state::ACPowerFlowState; bus_types=nothing)

Compute the 4-block power flow Jacobian at the current operating point.

Returns a NamedTuple with:

  • dp_dva: ∂P/∂θ (n × n) — J1
  • dp_dvm: ∂P/∂|V| (n × n) — J2
  • dq_dva: ∂Q/∂θ (n × n) — J3
  • dq_dvm: ∂Q/∂|V| (n × n) — J4

By default (bus_types=nothing), returns raw partial derivatives for ALL buses.

Pass bus_types::Vector{Int} (1=PQ, 2=PV, 3=slack) to apply Newton-Raphson bus-type column modifications matching PowerModels' calc_basic_jacobian_matrix:

  • PQ (1): Raw derivatives (no modification)
  • PV (2): θ columns unchanged; |V| columns zeroed with ∂Q_j/∂|V_j| = 1
  • Slack (3): All columns become unit vectors (e_j)
source

Utilities

PowerDiff.silenceFunction
silence()

Suppress all warning messages from PowerDiff for the rest of the session.

Warnings from other packages (PowerModels, JuMP, Ipopt, etc.) are not affected. To suppress PowerModels output, also call PowerModels.silence().

source

Abstract Types

PowerDiff.AbstractPowerNetworkType
AbstractPowerNetwork

Abstract base type for power network data containers.

Concrete subtypes:

  • DCNetwork: DC network with susceptance-weighted Laplacian
  • ACNetwork: AC network with complex admittance matrix
source
PowerDiff.AbstractPowerFlowStateType
AbstractPowerFlowState

Abstract base type for power flow solutions (DC or AC).

Concrete subtypes:

  • DCPowerFlowState: DC power flow solution (θr = Lr \ p_r)
  • ACPowerFlowState: AC power flow solution (complex voltages)
  • AbstractOPFSolution and its subtypes
source
PowerDiff.AbstractOPFSolutionType
AbstractOPFSolution <: AbstractPowerFlowState

Abstract type for optimal power flow solutions with dual variables.

Concrete subtypes:

  • DCOPFSolution: DC OPF with generation dispatch and duals
  • ACOPFSolution: AC OPF with voltages, generation, and duals
source
PowerDiff.AbstractOPFProblemType
AbstractOPFProblem

Abstract base type for OPF problem wrappers (JuMP models).

Concrete subtypes:

  • DCOPFProblem: DC OPF problem (B-θ formulation)
  • ACOPFProblem: AC OPF problem (polar formulation)
source

AcceleratedDCPowerFlows Extension

PowerDiff.to_apf_networkFunction
to_apf_network(net::DCNetwork) → APF.Network

Convert a DCNetwork to an AcceleratedDCPowerFlows network. Requires using AcceleratedDCPowerFlows.

source
PowerDiff.apf_ptdfFunction
apf_ptdf(net::DCNetwork; kwargs...) → APF.FullPTDF

Build an APF FullPTDF from a DCNetwork. Keyword arguments are forwarded to APF.full_ptdf. Requires using AcceleratedDCPowerFlows.

source
PowerDiff.apf_lodfFunction
apf_lodf(net::DCNetwork; kwargs...) → APF.FullLODF

Build an APF FullLODF from a DCNetwork. Keyword arguments are forwarded to APF.full_lodf. Requires using AcceleratedDCPowerFlows.

source
PowerDiff.compare_ptdfFunction
compare_ptdf(state::DCPowerFlowState; atol=1e-8) → NamedTuple

Cross-validate PowerDiff PTDF against AcceleratedDCPowerFlows PTDF. Returns (match, maxerr) where match is true if all entries agree within atol. Requires using AcceleratedDCPowerFlows.

source
PowerDiff.materialize_apf_ptdfFunction
materialize_apf_ptdf(Phi::APF.FullPTDF) → Matrix{Float64}

Materialize a dense PTDF matrix from an APF FullPTDF object. Requires using AcceleratedDCPowerFlows.

source