[de061d] | 1 | /*
|
---|
| 2 | * vmg - a versatile multigrid solver
|
---|
| 3 | * Copyright (C) 2012 Institute for Numerical Simulation, University of Bonn
|
---|
| 4 | *
|
---|
| 5 | * vmg is free software: you can redistribute it and/or modify
|
---|
| 6 | * it under the terms of the GNU General Public License as published by
|
---|
| 7 | * the Free Software Foundation, either version 3 of the License, or
|
---|
| 8 | * (at your option) any later version.
|
---|
| 9 | *
|
---|
| 10 | * vmg is distributed in the hope that it will be useful,
|
---|
| 11 | * but WITHOUT ANY WARRANTY; without even the implied warranty of
|
---|
| 12 | * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
---|
| 13 | * GNU General Public License for more details.
|
---|
| 14 | *
|
---|
| 15 | * You should have received a copy of the GNU General Public License
|
---|
| 16 | * along with this program. If not, see <http://www.gnu.org/licenses/>.
|
---|
| 17 | */
|
---|
| 18 |
|
---|
| 19 | /**
|
---|
| 20 | * @file givens.hpp
|
---|
| 21 | * @author Julian Iseringhausen <isering@ins.uni-bonn.de>
|
---|
| 22 | * @date Mon Apr 18 13:10:15 2011
|
---|
| 23 | *
|
---|
| 24 | * @brief Solves system of linear equations using Givens rotations.
|
---|
| 25 | * Compare to Meister, Numerik lineare Gleichungssysteme.
|
---|
| 26 | *
|
---|
| 27 | */
|
---|
| 28 |
|
---|
| 29 | #ifndef GIVENS_HPP_
|
---|
| 30 | #define GIVENS_HPP_
|
---|
| 31 |
|
---|
| 32 | #include <cfloat>
|
---|
| 33 |
|
---|
| 34 | #include "solver/solver.hpp"
|
---|
| 35 |
|
---|
| 36 | namespace VMG
|
---|
| 37 | {
|
---|
| 38 |
|
---|
| 39 | template<class T>
|
---|
| 40 | class Givens : public T
|
---|
| 41 | {
|
---|
| 42 | public:
|
---|
| 43 | Givens(bool register_ = true) :
|
---|
| 44 | T(register_)
|
---|
| 45 | {}
|
---|
| 46 |
|
---|
| 47 | Givens(int size, bool register_ = true) :
|
---|
| 48 | T(size, register_)
|
---|
| 49 | {}
|
---|
| 50 |
|
---|
| 51 | protected:
|
---|
| 52 | void Compute();
|
---|
| 53 | };
|
---|
| 54 |
|
---|
| 55 | template<class T>
|
---|
| 56 | void Givens<T>::Compute()
|
---|
| 57 | {
|
---|
| 58 | int n = this->Size();
|
---|
| 59 | vmg_float c,s,t;
|
---|
| 60 |
|
---|
| 61 | for (int i=0; i<n-1; i++)
|
---|
| 62 | for (int j=i+1; j<n; j++)
|
---|
| 63 | if (fabs(this->Mat(j,i)) > DBL_EPSILON) {
|
---|
| 64 |
|
---|
| 65 | t = 1.0 / sqrt(this->Mat(i,i)*this->Mat(i,i) + this->Mat(j,i)*this->Mat(j,i));
|
---|
| 66 | s = t * this->Mat(j,i);
|
---|
| 67 | c = t * this->Mat(i,i);
|
---|
| 68 |
|
---|
| 69 | for (int k=i; k<n; k++) {
|
---|
| 70 |
|
---|
| 71 | t = c * this->Mat(i,k) + s * this->Mat(j,k);
|
---|
| 72 |
|
---|
| 73 | if (k != i)
|
---|
| 74 | this->Mat(j,k) = c * this->Mat(j,k) - s * this->Mat(i,k);
|
---|
| 75 |
|
---|
| 76 | this->Mat(i,k) = t;
|
---|
| 77 |
|
---|
| 78 | }
|
---|
| 79 |
|
---|
| 80 | t = c * this->Rhs(i) + s * this->Rhs(j);
|
---|
| 81 |
|
---|
| 82 | this->Rhs(j) = c * this->Rhs(j) - s * this->Rhs(i);
|
---|
| 83 |
|
---|
| 84 | this->Rhs(i) = t;
|
---|
| 85 |
|
---|
| 86 | this->Mat(j,i) = 0.0;
|
---|
| 87 |
|
---|
| 88 | }
|
---|
| 89 |
|
---|
| 90 | for (int i=n-1; i>=0; i--) {
|
---|
| 91 |
|
---|
| 92 | for (int j=i+1; j<n; j++)
|
---|
| 93 | this->Rhs(i) -= this->Mat(i,j) * this->Sol(j);
|
---|
| 94 |
|
---|
| 95 | this->Sol(i) = this->Rhs(i) / this->Mat(i,i);
|
---|
| 96 |
|
---|
| 97 | }
|
---|
| 98 | }
|
---|
| 99 |
|
---|
| 100 | }
|
---|
| 101 |
|
---|
| 102 | #endif /* GIVENS_HPP_ */
|
---|