deltaFlow
Admittance.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 <complex>
27
28#include "Admittance.H"
29#include "Logger.H"
30#include "Data.H"
31
32Eigen::MatrixXcd computeAdmittanceMatrix(const BusData& busData, const BranchData& branchData) {
33 int nLine = branchData.From.size();
34 int N = std::max(branchData.From.maxCoeff(), branchData.To.maxCoeff()); // 1-based bus indexing
35
36 Eigen::MatrixXcd Ybus = Eigen::MatrixXcd::Zero(N, N); // Size N x N, 0-based indexing internally
37
38 for (int k = 0; k < nLine; ++k) {
39 int from = branchData.From(k) - 1; // Convert to 0-based
40 int to = branchData.To(k) - 1; // Convert to 0-based
41
42 std::complex<double> Z(branchData.R(k), branchData.X(k));
43 std::complex<double> Y = 1.0 / Z;
44
45 std::complex<double> B(0.0, 0.5 * branchData.B(k));
46
47 double a = branchData.tapRatio(k);
48 if (a == 0.0) {
49 a = 1.0;
50 }
51
52 // Off-diagonal
53 Ybus(from, to) -= Y / a;
54 Ybus(to, from) = Ybus(from, to); // Symmetric
55
56 // Diagonal
57 Ybus(from, from) += Y / (a * a) + B;
58 Ybus(to, to) += Y + B;
59 }
60
61 // Add shunt admittances
62 int nBuses = busData.ID.size(); // Should match Gs and Bs size
63
64 for (int n = 0; n < nBuses; ++n) {
65 int busIndex = busData.ID(n) - 1; // Convert to 0-based
66
67 if (busIndex < 0 || busIndex >= Ybus.rows()) {
68 LOG_ERROR("Warning: Bus ID {} out of bounds in Ybus", busIndex + 1);
69 continue;
70 }
71
72 std::complex<double> Y_shunt(busData.Gs(n), busData.Bs(n));
73 Ybus(busIndex, busIndex) += Y_shunt;
74 }
75
76 return Ybus;
77}
Eigen::MatrixXcd computeAdmittanceMatrix(const BusData &busData, const BranchData &branchData)
Computes the complex bus admittance matrix ($$ Y_{bus} $$).
Definition Admittance.C:32
Declaration of functions and classes for constructing the bus admittance matrix ($$ Y_{bus} $$) in po...
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_ERROR(msg,...)
Macro for logging an error-level message.
Definition Logger.H:88
Contains all relevant data for each transmission line or transformer branch.
Definition Data.H:83
Eigen::VectorXd tapRatio
Transformer tap ratio ($$ a $$)
Definition Data.H:90
Eigen::VectorXd B
Line susceptance ($$ B $$) [p.u.].
Definition Data.H:89
Eigen::VectorXd X
Reactance ($$ X $$) [p.u.].
Definition Data.H:87
Eigen::VectorXi From
From bus indices.
Definition Data.H:84
Eigen::VectorXd R
Resistance ($$ R $$) [p.u.].
Definition Data.H:86
Eigen::VectorXi To
To bus indices.
Definition Data.H:85
Contains all relevant data for each bus in the power system.
Definition Data.H:55
Eigen::VectorXi ID
Bus numbers.
Definition Data.H:56
Eigen::VectorXd Gs
Shunt conductance ($$ G_{sh} $$) [p.u.].
Definition Data.H:68
Eigen::VectorXd Bs
Shunt susceptance ($$ B_{sh} $$) [p.u.].
Definition Data.H:69