PSCF v1.1
LuSolver.h
1#ifndef PSCF_LU_SOLVER_H
2#define PSCF_LU_SOLVER_H
3
4/*
5* PSCF - Polymer Self-Consistent Field Theory
6*
7* Copyright 2016 - 2022, The Regents of the University of Minnesota
8* Distributed under the terms of the GNU General Public License.
9*/
10
11#include <util/containers/Array.h>
12#include <util/containers/Matrix.h>
13#include <gsl/gsl_matrix.h>
14#include <gsl/gsl_vector.h>
15#include <gsl/gsl_permutation.h>
16
17namespace Pscf
18{
19
20 using namespace Util;
21
31 {
32 public:
33
37 LuSolver();
38
42 ~LuSolver();
43
49 void allocate(int n);
50
56 void computeLU(const Matrix<double>& A);
57
65
71 void inverse (Matrix<double>& inv);
72
73 private:
74
76 gsl_vector b_;
77
79 gsl_vector x_;
80
82 gsl_matrix* luPtr_;
83
85 gsl_matrix* gMatInverse_;
86
88 gsl_permutation* permPtr_;
89
91 int signum_;
92
94 int n_;
95
96 };
97
98}
99#endif
Solve Ax=b by LU decomposition of A.
Definition: LuSolver.h:31
void inverse(Matrix< double > &inv)
Compute inverse of matrix A.
Definition: LuSolver.cpp:104
~LuSolver()
Destructor.
Definition: LuSolver.cpp:34
void allocate(int n)
Allocate memory.
Definition: LuSolver.cpp:46
void solve(Array< double > &b, Array< double > &x)
Solve Ax = b for known b to compute x.
Definition: LuSolver.cpp:83
LuSolver()
Constructor.
Definition: LuSolver.cpp:14
void computeLU(const Matrix< double > &A)
Compute the LU decomposition for later use.
Definition: LuSolver.cpp:63
Array container class template.
Definition: Array.h:34
Two-dimensional array container template (abstract).
Definition: Matrix.h:32
C++ namespace for polymer self-consistent field theory (PSCF).
Utility classes for scientific computation.
Definition: accumulators.mod:1