deltaFlow
Qlim.C File Reference

Reactive power limit enforcement implementation. More...

#include <cmath>
#include <limits>
#include <vector>
#include "Qlim.H"
#include "Logger.H"
#include "Data.H"
Include dependency graph for Qlim.C:

Go to the source code of this file.

Functions

bool checkQlimits (const Eigen::VectorXd &V, const Eigen::VectorXd &delta, Eigen::VectorXi &type_bus, const Eigen::MatrixXd &G, const Eigen::MatrixXd &B, BusData &busData, const std::vector< int > &pv_bus_id, int n_bus)
 Checks reactive power limits on PV buses after solver convergence.
 

Detailed Description

Reactive power limit enforcement implementation.

Definition in file Qlim.C.

Function Documentation

◆ checkQlimits()

bool checkQlimits ( const Eigen::VectorXd &  V,
const Eigen::VectorXd &  delta,
Eigen::VectorXi &  type_bus,
const Eigen::MatrixXd &  G,
const Eigen::MatrixXd &  B,
BusData busData,
const std::vector< int > &  pv_bus_id,
int  n_bus 
)

Checks reactive power limits on PV buses after solver convergence.

Computes the reactive power at each bus using converged voltages and angles, then checks if PV buses violate their $$ Q_{gmax} $$ or $$ Q_{gmin} $$ limits. Violating PV buses are converted to PQ type (type 3 in C++ convention).

Parameters
VConverged voltage magnitudes [p.u.].
deltaConverged voltage angles [rad].
type_bus(in/out) Bus type vector; PV buses that violate limits are set to PQ.
GConductance matrix.
BSusceptance matrix.
busDataBus data (for Ql, Qgmax, Qgmin).
pv_bus_id0-based indices of PV buses (before any conversion).
n_busTotal number of buses.
Returns
true if any Q-limit was hit (solver must be re-run), false otherwise.

Definition at line 34 of file Qlim.C.

43 {
44 // Q-limits: zero means no limit, using +/-inf
45 Eigen::VectorXd Qmax = busData.Qgmax;
46 Eigen::VectorXd Qmin = busData.Qgmin;
47
48 for (int i = 0; i < n_bus; ++i) {
49 if (Qmax(i) == 0.0) Qmax(i) = std::numeric_limits<double>::infinity();
50 if (Qmin(i) == 0.0) Qmin(i) = -std::numeric_limits<double>::infinity();
51 }
52
53 // Compute reactive power Q at each bus with converged V and delta
54 Eigen::VectorXd Q = Eigen::VectorXd::Zero(n_bus);
55 for (int i = 0; i < n_bus; ++i) {
56 for (int j = 0; j < n_bus; ++j) {
57 double dij = delta(i) - delta(j);
58 Q(i) += V(i) * V(j) * (G(i, j) * std::sin(dij) - B(i, j) * std::cos(dij));
59 }
60 }
61
62 // Qg = Q_calc + Ql (all in p.u.)
63 Eigen::VectorXd Qg = Q + busData.Ql;
64
65 // Check Q-limits only for PV buses
66 bool qlim_hit = false;
67
68 for (int idx : pv_bus_id) {
69 if (Qg(idx) > Qmax(idx)) {
70 type_bus(idx) = 3; // PV to PQ
71 busData.Qg(idx) = Qmax(idx); // Fix Q at limit for next solver run
72 qlim_hit = true;
73 LOG_DEBUG("Q-limit (max) hit at bus {} : Qg = {:.4f} > Qmax = {:.4f}", idx + 1, Qg(idx), Qmax(idx));
74 } else if (Qg(idx) < Qmin(idx)) {
75 type_bus(idx) = 3; // PV to PQ
76 busData.Qg(idx) = Qmin(idx); // Fix Q at limit for next solver run
77 qlim_hit = true;
78 LOG_DEBUG("Q-limit (min) hit at bus {} : Qg = {:.4f} < Qmin = {:.4f}", idx + 1, Qg(idx), Qmin(idx));
79 }
80 }
81
82 if (!qlim_hit) {
83 LOG_DEBUG("Power flow converged without hitting Q-limits.");
84 }
85
86 return qlim_hit;
87}
#define LOG_DEBUG(msg,...)
Macro for logging a debug-level message.
Definition Logger.H:85
Eigen::VectorXd Ql
Reactive power load [MVAr or p.u.].
Definition Data.H:65
Eigen::VectorXd Qgmax
Max reactive power generation [MVAr or p.u.].
Definition Data.H:66
Eigen::VectorXd Qgmin
Min reactive power generation [MVAr or p.u.].
Definition Data.H:67
Eigen::VectorXd Qg
Reactive power generation [MVAr or p.u.].
Definition Data.H:63

References LOG_DEBUG, BusData::Qg, BusData::Qgmax, BusData::Qgmin, and BusData::Ql.

Here is the caller graph for this function: