deltaFlow
PowerMismatch.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
28#include "PowerMismatch.H"
29
30Eigen::VectorXd powerMismatch(
31 const Eigen::VectorXd& Ps,
32 const Eigen::VectorXd& Qs,
33 const Eigen::MatrixXd& G,
34 const Eigen::MatrixXd& B,
35 const Eigen::VectorXd& V,
36 const Eigen::VectorXd& delta,
37 int n_bus,
38 const std::vector<int>& pq_bus_id,
39 Eigen::VectorXd& P,
40 Eigen::VectorXd& Q
41) {
42 P = Eigen::VectorXd::Zero(n_bus);
43 Q = Eigen::VectorXd::Zero(n_bus);
44
45 // Compute P(i) and Q(i) at each bus
46 for (int i = 0; i < n_bus; ++i) {
47 for (int j = 0; j < n_bus; ++j) {
48 double dij = delta(i) - delta(j);
49 P(i) += V(i) * V(j) * (G(i, j) * std::cos(dij) + B(i, j) * std::sin(dij));
50 Q(i) += V(i) * V(j) * (G(i, j) * std::sin(dij) - B(i, j) * std::cos(dij));
51 }
52 }
53
54 // delta_P for all non-slack buses (bus index 1..N-1, 0-based)
55 Eigen::VectorXd delta_P = Ps - P;
56 Eigen::VectorXd delta_Q = Qs - Q;
57
58 // Only keep delta_Q for PQ buses
59 int n_pq = static_cast<int>(pq_bus_id.size());
60 Eigen::VectorXd delta_Q_pq(n_pq);
61 for (int k = 0; k < n_pq; ++k) {
62 delta_Q_pq(k) = delta_Q(pq_bus_id[k]);
63 }
64
65 // Mismatch = [delta_P(2:end); delta_Q(pq)]
66 // In 0-based: delta_P(1..N-1)
67 Eigen::VectorXd mismatch(n_bus - 1 + n_pq);
68 mismatch.head(n_bus - 1) = delta_P.tail(n_bus - 1);
69 mismatch.tail(n_pq) = delta_Q_pq;
70
71 return mismatch;
72}
Eigen::VectorXd powerMismatch(const Eigen::VectorXd &Ps, const Eigen::VectorXd &Qs, const Eigen::MatrixXd &G, const Eigen::MatrixXd &B, const Eigen::VectorXd &V, const Eigen::VectorXd &delta, int n_bus, const std::vector< int > &pq_bus_id, Eigen::VectorXd &P, Eigen::VectorXd &Q)
Computes the power mismatch vector for use in Newton-Raphson iterations.
Power mismatch computation for Newton-Raphson power flow analysis.