PSCF v1.1
TridiagonalSolver.h
1#ifndef PSCF_TRIDIAGONAL_SOLVER_H
2#define PSCF_TRIDIAGONAL_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/DArray.h>
12
13#include <string>
14#include <iostream>
15
16
17namespace Pscf
18{
19
20 using namespace Util;
21
28 {
29 public:
30
35
40
46 void allocate(int n);
47
54 void computeLU(const DArray<double>& d, const DArray<double>& u);
55
63 void computeLU(const DArray<double>& d,
64 const DArray<double>& u,
65 const DArray<double>& l);
66
73 void multiply(const DArray<double>& b, DArray<double>& x);
74
81 void solve(const DArray<double>& b, DArray<double>& x);
82
83 private:
84
85 // Diagonal elements
87
88 // Upper off-diagonal elements (unmodified by computeLU)
90
91 // Lower off-diagonal elements (replaced by multipliers)
93
94 // Work space.
96
97 int n_;
98
99 // Apply Gauss elimination to private arrays d_, u_, l_.
100 void gaussElimination();
101
102 };
103
104}
105#endif
Solver for Ax=b with tridiagonal matrix A.
void multiply(const DArray< double > &b, DArray< double > &x)
Evaluate product Ab = x for known b to compute x.
void computeLU(const DArray< double > &d, const DArray< double > &u)
Compute LU decomposition of a symmetric tridiagonal matrix.
void solve(const DArray< double > &b, DArray< double > &x)
Solve Ax = b for known b to compute x.
void allocate(int n)
Allocate memory.
Dynamically allocatable contiguous array template.
Definition: DArray.h:32
C++ namespace for polymer self-consistent field theory (PSCF).
Utility classes for scientific computation.
Definition: accumulators.mod:1