deltaFlow
NewtonRaphson.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#include <limits>
28#include <iostream>
29#include <vector>
30
31#include "Logger.H"
32#include "NewtonRaphson.H"
33#include "Jacobian.H"
34#include "PowerMismatch.H"
35#include "Progress.H"
36#include "Data.H"
37#include "Utils.H"
38
40 const Eigen::MatrixXd& G,
41 const Eigen::MatrixXd& B,
42 const Eigen::VectorXd& Ps,
43 const Eigen::VectorXd& Qs,
44 Eigen::VectorXd& V,
45 Eigen::VectorXd& delta,
46 int n_bus,
47 int n_pq,
48 const std::vector<int>& pq_bus_id,
49 int maxIter,
50 double tolerance,
51 std::vector<std::pair<int, double>>* iterHistory
52) {
53 // Compute initial mismatch
54 Eigen::VectorXd P(n_bus), Q(n_bus);
55 Eigen::VectorXd mismatch = powerMismatch(Ps, Qs, G, B, V, delta, n_bus, pq_bus_id, P, Q);
56
57 double error = mismatch.cwiseAbs().maxCoeff();
58 int iter = 0;
59
60 if (iterHistory) {
61 iterHistory->clear();
62 iterHistory->emplace_back(0, error);
63 }
64
65 while (error >= tolerance) {
66 if (iter >= maxIter) {
67 printConvergenceStatus("Newton-Raphson", false, iter, maxIter, error, tolerance);
68 LOG_WARN("Newton-Raphson did not converge within {} iterations.", maxIter);
69 LOG_DEBUG("Final max mismatch was {:.6e}, tolerance is {:.6e}.", error, tolerance);
70 return false;
71 }
72 iter++;
73
74 // Build Jacobian
75 Eigen::MatrixXd J = computeJacobian(V, delta, n_bus, n_pq, pq_bus_id, G, B, P, Q);
76
77 // Solve J * correction = mismatch
78 Eigen::VectorXd correction = J.colPivHouseholderQr().solve(mismatch);
79
80 // Update delta for non-slack buses (indices 1..N-1)
81 for (int i = 1; i < n_bus; ++i) {
82 delta(i) += correction(i - 1);
83 }
84
85 // Update V for PQ buses only
86 for (int k = 0; k < n_pq; ++k) {
87 V(pq_bus_id[k]) += correction(n_bus - 1 + k);
88 }
89
90 // Recompute mismatch
91 mismatch = powerMismatch(Ps, Qs, G, B, V, delta, n_bus, pq_bus_id, P, Q);
92
93 error = mismatch.cwiseAbs().maxCoeff();
94 if (iterHistory) iterHistory->emplace_back(iter, error);
95 printIterationProgress("Newton-Raphson", iter, maxIter, error, tolerance);
96 LOG_DEBUG("NR iteration {}: max mismatch = {:.16e}", iter, error);
97 }
98
99 printConvergenceStatus("Newton-Raphson", true, iter, maxIter, error, tolerance);
100 LOG_DEBUG("Newton-Raphson converged in {} iterations with max mismatch {:.6e}", iter, error);
101 return true;
102}
@ NewtonRaphson
Newton-Raphson iterative method.
Data structures and utility functions for reading and displaying bus and branch data in power system ...
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.
Definition Jacobian.C:30
Jacobian matrix computation for Newton-Raphson power flow analysis.
Logger utilities for deltaFlow, providing logging macros and a singleton Logger class.
#define LOG_WARN(msg,...)
Macro for logging a warning-level message.
Definition Logger.H:87
#define LOG_DEBUG(msg,...)
Macro for logging a debug-level message.
Definition Logger.H:85
Declaration of the Newton-Raphson load flow solver for power system analysis.
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.
Progress bar for iterative solvers.
void printIterationProgress(const std::string &solver, int iter, int maxIter, double error, double tol, int barWidth=50)
Print an iteration progress line for a solver.
Definition Progress.H:50
void printConvergenceStatus(const std::string &solver, bool converged, int iter, int maxIter, double error, double tol, int barWidth=50)
Print the final convergence status line.
Definition Progress.H:100
Utility functions and helpers for deltaFlow.