deltaFlow
Jacobian.H File Reference

Jacobian matrix computation for Newton-Raphson power flow analysis. More...

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

Go to the source code of this file.

Functions

Eigen::MatrixXd computeJacobian (const Eigen::VectorXd &V, const Eigen::VectorXd &delta, int n_bus, int n_pq, const std::vector< int > &pq_bus_id, const Eigen::MatrixXd &G, const Eigen::MatrixXd &B, const Eigen::VectorXd &P, const Eigen::VectorXd &Q)
 Computes the Jacobian matrix for the Newton-Raphson power flow solver.
 

Detailed Description

Jacobian matrix computation for Newton-Raphson power flow analysis.

Constructs the Jacobian matrix with block structure:

$$ J = \begin{bmatrix} J_{11} & J_{12} \ J_{21} & J_{22} \end{bmatrix} $$

where:

  • $$ J_{11} $$: $$ \frac{\partial P}{\partial \delta} $$ for non-slack buses
  • $$ J_{12} $$: $$ \frac{\partial P}{\partial |V|} $$ for non-slack rows, PQ columns
  • $$ J_{21} $$: $$ \frac{\partial Q}{\partial \delta} $$ for PQ rows, non-slack columns
  • $$ J_{22} $$: $$ \frac{\partial Q}{\partial |V|} $$ for PQ rows and columns

Follows the formulation by B. Sereeter, C. Vuik, and C. Witteveen (REPORT 17-07).

Definition in file Jacobian.H.

Function Documentation

◆ computeJacobian()

Eigen::MatrixXd computeJacobian ( const Eigen::VectorXd &  V,
const Eigen::VectorXd &  delta,
int  n_bus,
int  n_pq,
const std::vector< int > &  pq_bus_id,
const Eigen::MatrixXd &  G,
const Eigen::MatrixXd &  B,
const Eigen::VectorXd &  P,
const Eigen::VectorXd &  Q 
)

Computes the Jacobian matrix for the Newton-Raphson power flow solver.

Parameters
VVoltage magnitudes at each bus [p.u.].
deltaVoltage angles at each bus [rad].
n_busTotal number of buses.
n_pqNumber of PQ buses.
pq_bus_id0-based indices of PQ buses.
GConductance matrix (real part of $$ Y_{bus} $$).
BSusceptance matrix (imaginary part of $$ Y_{bus} $$).
PComputed active power at each bus (from power mismatch) [p.u.].
QComputed reactive power at each bus (from power mismatch) [p.u.].
Returns
The assembled Jacobian matrix.

Definition at line 30 of file Jacobian.C.

40 {
41 // J11: (n_bus-1) x (n_bus-1) - dP/ddelta, non-slack buses
42 Eigen::MatrixXd J11 = Eigen::MatrixXd::Zero(n_bus - 1, n_bus - 1);
43
44 for (int i = 1; i < n_bus; ++i) {
45 for (int k = 1; k < n_bus; ++k) {
46 if (i == k) {
47 J11(i - 1, k - 1) = -Q(i) - V(i) * V(i) * B(i, i);
48 } else {
49 double dik = delta(i) - delta(k);
50 J11(i - 1, k - 1) = V(i) * V(k) * (G(i, k) * std::sin(dik)
51 - B(i, k) * std::cos(dik));
52 }
53 }
54 }
55
56 // J21: n_pq x (n_bus-1) - dQ/ddelta, PQ rows, non-slack columns
57 Eigen::MatrixXd J21 = Eigen::MatrixXd::Zero(n_pq, n_bus - 1);
58
59 for (int i = 0; i < n_pq; ++i) {
60 int j = pq_bus_id[i];
61 for (int k = 1; k < n_bus; ++k) {
62 if (j == k) {
63 J21(i, k - 1) = P(j) - V(j) * V(j) * G(j, j);
64 } else {
65 double djk = delta(j) - delta(k);
66 J21(i, k - 1) = -V(j) * V(k) * (G(j, k) * std::cos(djk)
67 + B(j, k) * std::sin(djk));
68 }
69 }
70 }
71
72 // J12: (n_bus-1) x n_pq - dP/dV, non-slack rows, PQ columns
73 Eigen::MatrixXd J12 = Eigen::MatrixXd::Zero(n_bus - 1, n_pq);
74
75 for (int i = 1; i < n_bus; ++i) {
76 for (int k = 0; k < n_pq; ++k) {
77 int j = pq_bus_id[k];
78 if (i == j) {
79 J12(i - 1, k) = P(j) / V(j) + V(j) * G(j, j);
80 } else {
81 double dij = delta(i) - delta(j);
82 J12(i - 1, k) = V(i) * (G(i, j) * std::cos(dij)
83 + B(i, j) * std::sin(dij));
84 }
85 }
86 }
87
88 // J22: n_pq x n_pq - dQ/dV, PQ rows and columns
89 Eigen::MatrixXd J22 = Eigen::MatrixXd::Zero(n_pq, n_pq);
90
91 for (int i = 0; i < n_pq; ++i) {
92 int j = pq_bus_id[i];
93 for (int k = 0; k < n_pq; ++k) {
94 int l = pq_bus_id[k];
95 if (j == l) {
96 J22(i, k) = Q(j) / V(j) - V(j) * B(j, j);
97 } else {
98 double djl = delta(j) - delta(l);
99 J22(i, k) = V(j) * (G(j, l) * std::sin(djl)
100 - B(j, l) * std::cos(djl));
101 }
102 }
103 }
104
105 // Assemble full Jacobian: [J11 J12; J21 J22]
106 int dim = (n_bus - 1) + n_pq;
107 Eigen::MatrixXd J(dim, dim);
108 J.topLeftCorner(n_bus - 1, n_bus - 1) = J11;
109 J.topRightCorner(n_bus - 1, n_pq) = J12;
110 J.bottomLeftCorner(n_pq, n_bus - 1) = J21;
111 J.bottomRightCorner(n_pq, n_pq) = J22;
112
113 return J;
114}
Here is the caller graph for this function: