deltaFlow
PowerMismatch.C File Reference

Active and reactive power mismatch computation implementation. More...

#include <cmath>
#include "PowerMismatch.H"
Include dependency graph for PowerMismatch.C:

Go to the source code of this file.

Functions

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.
 

Detailed Description

Active and reactive power mismatch computation implementation.

Definition in file PowerMismatch.C.

Function Documentation

◆ powerMismatch()

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.

The mismatch vector is composed of:

  • $$ \Delta P $$ for all non-slack buses (buses 1..N-1, 0-indexed)
  • $$ \Delta Q $$ for PQ buses only
Parameters
PsScheduled active power injection (Pg - Pl) at each bus [p.u.].
QsScheduled reactive power injection (Qg - Ql) at each bus [p.u.].
GConductance matrix (real part of $$ Y_{bus} $$).
BSusceptance matrix (imaginary part of $$ Y_{bus} $$).
VVoltage magnitudes at each bus [p.u.].
deltaVoltage angles at each bus [rad].
n_busTotal number of buses.
pq_bus_id0-based indices of PQ buses.
P(output) Computed active power at each bus [p.u.].
Q(output) Computed reactive power at each bus [p.u.].
Returns
Mismatch vector $$ [\Delta P_{2..N}; \Delta Q_{PQ}] $$.

Definition at line 30 of file PowerMismatch.C.

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}
Here is the caller graph for this function: