xlu.h
Go to the documentation of this file.
1 #pragma once
2 
3 #include "xmatrix.h"
4 
5 namespace nmx::gsl {
6 
12 {
13 public:
16 
17 private:
18  //linke Seite des Gleichungssystems
19  Matrix _left;
20  //rechte Seite des Gleichungssystems
21  Vector _right;
22  // Lösungsvektor
23  Vector _solution;
24 
25 public:
26  //Reihen- und Zeilenvertauschung (gsl)
27  gsl_permutation *_permutation;
28  // (-1)^n Anzahl der Vertauschungen (gsl)
30 
31 public:
35  inline LU_Decomposition() {}
36 
42  inline LU_Decomposition(const Matrix &m, const Vector &v)
43  : _left{ m }
44  , _right{ v }
45  , _solution(_right.size()) {}
46 
50  inline ~LU_Decomposition() {
51  if (_permutation != nullptr) {
52  //gsl-Struktur wird freigegeben
53  gsl_permutation_free(_permutation);
54  }
55  }
56 
61  inline void set_left(const Matrix &m) { _left = m; }
62 
67  inline void set_right(const Vector &v) { _right = v; }
68 
73  inline const Matrix &get_left() const { return _left; }
74 
80  const Vector &get_right() const { return _right; }
81 
87  inline double operator()(size_t idx) const { //
88  return _solution[idx];
89  }
90 
97  inline double operator()(size_t i, size_t j) const { //
98  return _left(i, j);
99  }
100 
105  inline void decompose() {
106  _permutation = gsl_permutation_alloc(std::min(_left.rows(), //
107  _left.columns()));
108  gsl_linalg_LU_decomp(_left.gsl(),
109  _permutation, //
111  }
112 
116  inline void solve() {
117  //zerlege Koeffizientenmatrix in obere und untere Dreiecksmatrix
118  decompose();
119  //Lösungsvektor muss die gleiche Länge haben wie die
120  //rechte Seite
121  if (_solution.size() != _right.size()) {
122  _solution = Vector(_right.size());
123  }
124  //löse das Gleichungssystem
125  gsl_linalg_LU_solve(_left.gsl(), //
126  _permutation,
127  _right.gsl(),
128  _solution.gsl());
129  }
130 
137  inline void solve(const Matrix &m, const Vector &v) {
138  set_left(m);
139  set_right(v);
140  solve();
141  }
142 
147  inline const Vector &solution() const { return _solution; }
148 
153  inline Matrix lower_matrix() const {
154  //reserviere Speicher fur das Ergebnis
155  Matrix mout(_left.rows(), _left.columns());
156 
157  for (size_t i = 0; i < _left.rows(); i++) {
158  for (size_t j = 0; j < _left.columns(); j++) {
159  if (i < j) {
160  mout(i, j) = 0;
161  } else if (i == j) {
162  mout(i, j) = 1.0; //diagonale
163  } else
164  mout(i, j) = _left(i, j);
165  }
166  }
167 
168  return mout;
169  }
170 
175  inline Matrix upper_matrix() const {
176  //reserviere Speicher fur das Ergebnis
177  Matrix mout(_left.rows(), _left.columns());
178 
179  for (size_t i = 0; i < _left.rows(); i++) {
180  for (size_t j = 0; j < _left.columns(); j++) {
181  if (i > j) {
182  mout(i, j) = 0;
183  } else
184  mout(i, j) = _left(i, j);
185  }
186  }
187 
188  return mout;
189  }
190 
196  //Anzahl der Elemente der Permutation
197  size_t psize = static_cast<size_t>(_permutation->size);
198 
199  Matrix pMatrix(_left.rows(), _left.columns());
200  for (size_t idx = 0; idx < psize; idx++) {
201  //für jede Reihe gebe die Position des Elements,
202  //welches 1 ist
203  size_t pIdx = gsl_permutation_get(_permutation, idx);
204  pMatrix(idx, pIdx) = 1;
205  }
206  return pMatrix;
207  }
208 
214  inline static double det(const Matrix &m) {
215  LU_Decomposition lu;
216  lu._left = m;
217  //LU-Zerlegung der Matrix
218  lu.decompose();
219  //Berechnung der Determinante
220  return gsl_linalg_LU_det(lu._left.gsl(), lu._permutation_sign);
221  }
222 
228  inline static Matrix inverse(const Matrix &m) {
229  //reserviere Speicher für die Inverse Matrix
230  Matrix mInverse(m.rows(), m.columns());
231  LU_Decomposition lu;
232  lu._left = m;
233  //LU-Zerlegung der Matrix
234  lu.decompose();
235  //Berechnung der Inversen Matrix
236  gsl_linalg_LU_invert(lu.get_left().gsl(), //
237  lu._permutation,
238  mInverse.gsl());
239  return mInverse;
240  }
241 
247  inline Vector permutation_vector() const {
248  //erzeuge Vektor für die Elemente der Permutation
249  gsl::Vector vec(_permutation->size);
250 
251  for (size_t idx = 0; idx < vec.size(); idx++) {
252  vec[idx] = gsl_permutation_get(_permutation, idx);
253  }
254  return vec;
255  }
256 };
257 
258 } // namespace nmx::gsl
const gsl_vector * gsl() const
gsl ermöglicht direkten Zugriff auf gsl-Funktionen
Definition: xvector.h:192
The Matrix class Schnittstelle zur gsl Bibliothek.
Definition: xmatrix.h:55
Vector permutation_vector() const
permutation_vector Permutationsvektor wird aus Permutationstruktur generiert
Definition: xlu.h:247
Matrix upper_matrix() const
upper_matrix
Definition: xlu.h:175
void decompose()
decompose Zerlegung der Koeffizientenmatrix in obere und untere Dreiecksmatrix
Definition: xlu.h:105
void set_right(const Vector &v)
set_right
Definition: xlu.h:67
size_t size() const
size
Definition: xvector1.h:25
void set_left(const Matrix &m)
set_left
Definition: xlu.h:61
Matrix lower_matrix() const
lower_matrix
Definition: xlu.h:153
gsl_permutation * _permutation
Definition: xlu.h:27
void solve(const Matrix &m, const Vector &v)
solve Berechnung des Lösungsvektors eines Gleichungssystems
Definition: xlu.h:137
LU::Vector Vector
Definition: x038.h:10
void solve()
solve Berechnung des Lösungsvektors
Definition: xlu.h:116
const Matrix & get_left() const
get_left
Definition: xlu.h:73
size_t columns() const
columns
Definition: xmatrix0.h:30
LU_Decomposition()
LUDecomposition erzeuge eine leere Instanz.
Definition: xlu.h:35
LU::Matrix Matrix
Definition: x038.h:9
size_t rows() const
rows
Definition: xmatrix0.h:24
double operator()(size_t i, size_t j) const
operator () Zugriff auf die Koeffizientenmatrix
Definition: xlu.h:97
The LU_Decomposition class Lösung eines linearen Gleichungssystems mittels der LU-Zerlegung.
Definition: xlu.h:11
The Vector class Klasse für gsl-Vektoren zussamengesetzt aus Komponenten.
Definition: xvector.h:76
LU_Decomposition(const Matrix &m, const Vector &v)
LUDecomposition Konstruktor.
Definition: xlu.h:42
Matrix permutation_matrix()
permutation_matrix
Definition: xlu.h:195
const Vector & solution() const
solution
Definition: xlu.h:147
double operator()(size_t idx) const
operator () Zugriff auf die Elemente des Lösungsvektors
Definition: xlu.h:87
const gsl_matrix * gsl() const
gsl_vec ermöglicht direkten Zugriff auf gsl Funktionen
Definition: xmatrix.h:161
const Vector & get_right() const
get_right
Definition: xlu.h:80
static Matrix inverse(const Matrix &m)
inverse Berechnung der Inversen einer Matrix
Definition: xlu.h:228
gsl::Vector Vector
Definition: xlu.h:15
static double det(const Matrix &m)
det Berechnung der Determinante
Definition: xlu.h:214