deltaFlow
Admittance.H File Reference

Declaration of functions and classes for constructing the bus admittance matrix ($$ Y_{bus} $$) in power system analysis. More...

#include <Eigen/Dense>
Include dependency graph for Admittance.H:
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Functions

Eigen::MatrixXcd computeAdmittanceMatrix (const BusData &busData, const BranchData &branchData)
 Computes the complex bus admittance matrix ($$ Y_{bus} $$).
 

Detailed Description

Declaration of functions and classes for constructing the bus admittance matrix ($$ Y_{bus} $$) in power system analysis.

This header provides the interface for computing the admittance matrix using bus and branch data. The admittance matrix $$ Y_{bus} $$ is a fundamental element in load flow studies, relating bus currents $$ I $$ and bus voltages $$ V $$ via the equation:

$$ I = Y_{bus} V $$

where:

  • $$ I $$ is the vector of injected currents,
  • $$ Y_{bus} $$ is the complex bus admittance matrix,
  • $$ V $$ is the vector of bus voltages.

Definition in file Admittance.H.

Function Documentation

◆ computeAdmittanceMatrix()

Eigen::MatrixXcd computeAdmittanceMatrix ( const BusData busData,
const BranchData branchData 
)

Computes the complex bus admittance matrix ($$ Y_{bus} $$).

This function generates the bus admittance matrix, which is essential in solving power flow equations in electrical networks. The admittance matrix is computed using the provided bus and branch data.

Parameters
busDataData representing the buses in the network.
branchDataData representing the branches (lines/transformers) in the network.
Returns
Eigen::MatrixXcd The computed $$ Y_{bus} $$ matrix as a complex matrix.

Definition at line 32 of file Admittance.C.

32 {
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}
#define LOG_ERROR(msg,...)
Macro for logging an error-level message.
Definition Logger.H:88
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
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

References BranchData::B, BusData::Bs, BranchData::From, BusData::Gs, BusData::ID, LOG_ERROR, BranchData::R, BranchData::tapRatio, BranchData::To, and BranchData::X.

Here is the caller graph for this function: