TBCI Numerical high perf. C++ Library 2.8.0
diluprecond.h
Go to the documentation of this file.
1
11
12#ifndef TBCI_SOLVER_DILUPRECOND_H
13#define TBCI_SOLVER_DILUPRECOND_H
14
15#include "tbci/solver/precond.h"
16
18
19// ********************************************************************************
20// ... For general BdMatrices
21// ********************************************************************************
22
23
24template <typename T> class BdMatrix;
25
39template <typename T>
40class DILU_BdMatrixPreconditioner : public Preconditioner_Sig<T, BdMatrix<T> >
41{
42 public:
47
48 // update function
49 void update (const BdMatrix<T> &A);
50
51 // solve functions
53
54 inline TVector<T> solve (const Vector<T> &v) const
55 {
56 TVector<T> tv(v);
57 return (solve(tv));
58 }
59
60 /* template <typename T> */
61 inline TVector<T> transSolve (const Vector<T> &v) const
62 { return solve(v); }
63 /* template <typename T> */
64 inline TVector<T> transSolve ( TVector<T> tv) const
65 { return solve(tv); }
66
67 private:
68 BVector<T> pivots;
69 unsigned int dimension;
70};
71
72// See Templates Book, Fig. 3.3
73template <typename T>
75{
76 // copy matrix variables to private vars
77 dimension = A.dim;
78
79 // calculate Pivots
80 if (pivots.size() != dimension)
81 pivots.resize(dimension);
82 BCHK(pivots.size() == 0, VecErr, "Null ptr in DILUBdMPrecCond::update", 0, );
83 const unsigned ds = A.diagconf.size();
84 unsigned i;
85 for (i = 0; i < dimension; ++i)
86 pivots.set(i) = A.diag[i];
87 for (i = 0; i < dimension; ++i) {
88 REGISTER T tmp = pivots.get(i);
89 if (UNLIKELY(tmp == T(0.0)))
90 tmp = T(1.0);
91 else
92 tmp = T(1.0)/tmp;
93 pivots.set(i) = tmp;
94 for (unsigned j = 0; j < ds; j++) {
95 const int df = A.diagconf.get(j);
96 if (LIKELY(df < (signed)dimension - (signed)i))
97 pivots.set(i+df) -= A.adiag.get(df)[i]*tmp*A.bdiag.get(df)[i];
98 }
99 }
100}
101
102template <typename T>
104{
105
106 // check dimensions
107 BCHKNR(dimension != y.size(), VecErr, "DILU_BdMPreCond::solve: Dimension mismatch", y.size());
108
109 // check trivial solution
110 // Matrix consists of main diagonal only
111 for (unsigned int i = 0; i < dimension; ++i)
112 y.setval(i) *= pivots.get(i);
113 return y;
114}
115
117
118#endif /* TBCI_SOLVER_DILUPRECOND_H */
const Vector< T > const Vector< T > & x
Definition LM_fit.h:97
int i
Definition LM_fit.h:71
doublereal y
Definition TOMS_707.C:27
#define BCHK(cond, exc, txt, ind, rtval)
Definition basics.h:575
#define LIKELY(expr)
branch prediction note that we sometimes on purpose mark the unlikely possibility likely and vice ver...
Definition basics.h:100
#define NAMESPACE_END
Definition basics.h:323
#define NAMESPACE_TBCI
Definition basics.h:317
#define BCHKNR(cond, exc, txt, ind)
Definition basics.h:586
#define UNLIKELY(expr)
Definition basics.h:101
#define REGISTER
Definition basics.h:108
#define T
Definition bdmatlib.cc:20
provides basic Vector functionality but arithmetic operators (+=, - , *, /...).
Definition bvector.h:68
unsigned long size() const HOT
Definition bvector.h:144
tbci_traits< T >::const_refval_type get(const unsigned long idx) const HOT
Definition bvector.h:132
The class BdMatrix is an implementation to store and do operations on sparse Matrices with a band str...
BVector< T * > adiag
BVector< T * > bdiag
pointers to upper and lower diagonals
unsigned int dim
BVector< unsigned > diagconf
configuration
TVector< T > transSolve(TVector< T > tv) const
Definition diluprecond.h:64
DILU_BdMatrixPreconditioner(const BdMatrix< T > &A)
Definition diluprecond.h:44
void update(const BdMatrix< T > &A)
Definition diluprecond.h:74
TVector< T > transSolve(const Vector< T > &v) const
Definition diluprecond.h:61
TVector< T > solve(TVector< T > x) const
TVector< T > solve(const Vector< T > &v) const
Definition diluprecond.h:54
Temporary Base Class Idiom: Class TVector is used for temporary variables.
Definition vector.h:73
exception class
Definition bvector.h:32