deltaFlow
Qlim.C
Go to the documentation of this file.
1/*
2 * Copyright (c) 2024 Saud Zahir
3 *
4 * This file is part of deltaFlow.
5 *
6 * deltaFlow is free software; you can redistribute it and/or
7 * modify it under the terms of the GNU General Public
8 * License as published by the Free Software Foundation; either
9 * version 3 of the License, or (at your option) any later version.
10 *
11 * deltaFlow is distributed in the hope that it will be useful,
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
14 * General Public License for more details.
15 *
16 * You should have received a copy of the GNU General Public
17 * License along with deltaFlow. If not, see
18 * <https://www.gnu.org/licenses/>.
19 */
20
26#include <cmath>
27#include <limits>
28#include <vector>
29
30#include "Qlim.H"
31#include "Logger.H"
32#include "Data.H"
33
35 const Eigen::VectorXd& V,
36 const Eigen::VectorXd& delta,
37 Eigen::VectorXi& type_bus,
38 const Eigen::MatrixXd& G,
39 const Eigen::MatrixXd& B,
40 BusData& busData,
41 const std::vector<int>& pv_bus_id,
42 int n_bus
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}
Data structures and utility functions for reading and displaying bus and branch data in power system ...
Logger utilities for deltaFlow, providing logging macros and a singleton Logger class.
#define LOG_DEBUG(msg,...)
Macro for logging a debug-level message.
Definition Logger.H:85
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.
Definition Qlim.C:34
Reactive power limit (Q-limit) checking for PV buses.
Contains all relevant data for each bus in the power system.
Definition Data.H:55
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