G+Smo  25.01.0
Geometry + Simulation Modules
 
Loading...
Searching...
No Matches
gsSimplePreconditioners.h
Go to the documentation of this file.
1
13#pragma once
14
17
18namespace gismo
19{
20
21namespace internal
22{
23template<typename T>
24void gaussSeidelSweep(const gsSparseMatrix<T> & A, gsMatrix<T>& x, const gsMatrix<T>& f);
25template<typename T>
26void reverseGaussSeidelSweep(const gsSparseMatrix<T> & A, gsMatrix<T>& x, const gsMatrix<T>& f);
27} // namespace internal
28
32template <typename MatrixType>
33class gsRichardsonOp GISMO_FINAL : public gsPreconditionerOp<typename MatrixType::Scalar>
34{
35 typedef memory::shared_ptr<MatrixType> MatrixPtr;
36 typedef typename MatrixType::Nested NestedMatrix;
37
38public:
40 typedef typename MatrixType::Scalar T;
41
43 typedef memory::shared_ptr< gsRichardsonOp > Ptr;
44
46 typedef memory::unique_ptr< gsRichardsonOp > uPtr;
47
50
52 explicit gsRichardsonOp(const MatrixType& mat, T tau = 1)
53 : m_mat(), m_expr(mat.derived()), m_tau(tau) {}
54
56 explicit gsRichardsonOp(const MatrixPtr& mat, T tau = 1)
57 : m_mat(mat), m_expr(m_mat->derived()), m_tau(tau) { }
58
59 static uPtr make(const MatrixType& mat, T tau = 1)
60 { return memory::make_unique( new gsRichardsonOp(mat, tau) ); }
61
62 static uPtr make(const MatrixPtr& mat, T tau = 1)
63 { return memory::make_unique( new gsRichardsonOp(mat, tau) ); }
64
65 void step(const gsMatrix<T> & rhs, gsMatrix<T> & x) const
66 {
67 GISMO_ASSERT( m_expr.rows() == rhs.rows() && m_expr.cols() == m_expr.rows(),
68 "Dimensions do not match.");
69
70 x += m_tau * ( rhs - m_expr * x );
71 }
72
73 // We use our own apply implementation as we can save one multiplication. This is important if the number
74 // of sweeps is 1: Then we can save *all* multiplications.
75 void apply(const gsMatrix<T> & input, gsMatrix<T> & x) const
76 {
77 GISMO_ASSERT( m_expr.rows() == input.rows() && m_expr.cols() == m_expr.rows(),
78 "Dimensions do not match.");
79
80 // For the first sweep, we do not need to multiply with the matrix
81 x.noalias() = m_tau * input;
82
83 for (index_t k = 1; k < m_num_of_sweeps; ++k)
84 x += m_tau * ( input - m_expr * x );
85 }
86
87 index_t rows() const {return m_expr.rows();}
88 index_t cols() const {return m_expr.cols();}
89
91 void setDamping(const T tau) { m_tau = tau; }
92
94 void getDamping() { return m_tau; }
95
98 {
100 opt.addReal( "Damping", "Damping parameter of the Richardson iteration", 1 );
101 return opt;
102 }
103
105 virtual void setOptions(const gsOptionList & opt)
106 {
107 Base::setOptions(opt);
108 m_tau = opt.askReal( "Damping", m_tau );
109 }
110
112 NestedMatrix matrix() const { return m_expr; }
113
115 MatrixPtr matrixPtr() const {
116 GISMO_ENSURE( m_mat, "A shared pointer is only available if it was provided to gsRichardsonOp." );
117 return m_mat;
118 }
119
120 typename gsLinearOperator<T>::Ptr underlyingOp() const { return makeMatrixOp(m_mat); }
121
122private:
123 const MatrixPtr m_mat;
124 NestedMatrix m_expr;
125
126 using Base::m_num_of_sweeps;
127 T m_tau;
128};
129
132template <class Derived>
133typename gsRichardsonOp<Derived>::uPtr makeRichardsonOp(const gsEigen::EigenBase<Derived>& mat, typename Derived::Scalar tau = 1)
134{ return gsRichardsonOp<Derived>::make(mat.derived(), tau); }
135
138template <class Derived>
139typename gsRichardsonOp<Derived>::uPtr makeRichardsonOp(const memory::shared_ptr<Derived>& mat, typename Derived::Scalar tau = 1)
140{ return gsRichardsonOp<Derived>::make(mat, tau); }
141
145template <typename MatrixType>
146class gsJacobiOp GISMO_FINAL : public gsPreconditionerOp<typename MatrixType::Scalar>
147{
148 typedef memory::shared_ptr<MatrixType> MatrixPtr;
149 typedef typename MatrixType::Nested NestedMatrix;
150
151public:
153 typedef typename MatrixType::Scalar T;
154
156 typedef memory::shared_ptr< gsJacobiOp > Ptr;
157
159 typedef memory::unique_ptr< gsJacobiOp > uPtr;
160
163
165 explicit gsJacobiOp(const MatrixType& mat, T tau = 1)
166 : m_mat(), m_expr(mat.derived()), m_tau(tau) {}
167
169 explicit gsJacobiOp(const MatrixPtr& mat, T tau = 1)
170 : m_mat(mat), m_expr(m_mat->derived()), m_tau(tau) { }
171
172 static uPtr make(const MatrixType& mat, T tau = 1)
173 { return memory::make_unique( new gsJacobiOp(mat, tau) ); }
174
175 static uPtr make(const MatrixPtr& mat, T tau = 1)
176 { return memory::make_unique( new gsJacobiOp(mat, tau) ); }
177
178 void step(const gsMatrix<T> & rhs, gsMatrix<T> & x) const
179 {
180 GISMO_ASSERT( m_expr.rows() == rhs.rows() && m_expr.cols() == m_expr.rows(),
181 "Dimensions do not match.");
182
183 GISMO_ASSERT( rhs.cols() == 1, "This operator is only implemented for a single right-hand side." );
184
185 x.array() += m_tau * ( rhs - m_expr * x ).array() / m_expr.diagonal().array();
186 }
187
188 // We use our own apply implementation as we can save one multiplication. This is important if the number
189 // of sweeps is 1: Then we can save *all* multiplications.
190 void apply(const gsMatrix<T> & input, gsMatrix<T> & x) const
191 {
192 GISMO_ASSERT( m_expr.rows() == input.rows() && m_expr.cols() == m_expr.rows(),
193 "Dimensions do not match.");
194
195 GISMO_ASSERT( input.cols() == 1, "This operator is only implemented for a single right-hand side." );
196
197 // For the first sweep, we do not need to multiply with the matrix
198 x.array() = m_tau * input.array() / m_expr.diagonal().array();
199
200 for (index_t k = 1; k < m_num_of_sweeps; ++k)
201 x.array() += m_tau * ( input - m_expr * x ).array() / m_expr.diagonal().array();
202 }
203
204 index_t rows() const {return m_expr.rows();}
205 index_t cols() const {return m_expr.cols();}
206
208 void setDamping(const T tau) { m_tau = tau; }
209
211 void getDamping() { return m_tau; }
212
215 {
217 opt.addReal( "Damping", "Damping parameter of the Jacobi iteration", 1 );
218 return opt;
219 }
220
222 virtual void setOptions(const gsOptionList & opt)
223 {
224 Base::setOptions(opt);
225 m_tau = opt.askReal( "Damping", m_tau );
226 }
227
229 NestedMatrix matrix() const { return m_expr; }
230
232 MatrixPtr matrixPtr() const {
233 GISMO_ENSURE( m_mat, "A shared pointer is only available if it was provided to gsJacobiOp." );
234 return m_mat;
235 }
236
237 typename gsLinearOperator<T>::Ptr underlyingOp() const { return makeMatrixOp(m_mat); }
238
239private:
240 const MatrixPtr m_mat;
241 NestedMatrix m_expr;
242 using Base::m_num_of_sweeps;
243 T m_tau;
244};
245
246
249template <class Derived>
250typename gsJacobiOp<Derived>::uPtr makeJacobiOp(const gsEigen::EigenBase<Derived>& mat, typename Derived::Scalar tau = 1)
251{ return gsJacobiOp<Derived>::make(mat.derived(), tau); }
252
255template <class Derived>
256typename gsJacobiOp<Derived>::uPtr makeJacobiOp(const memory::shared_ptr<Derived>& mat, typename Derived::Scalar tau = 1)
257{ return gsJacobiOp<Derived>::make(mat, tau); }
258
259
260struct gsGaussSeidel
261{
265 {
268 symmetric = 2
269 };
270};
271
277template <typename MatrixType, gsGaussSeidel::ordering ordering = gsGaussSeidel::forward>
278class gsGaussSeidelOp GISMO_FINAL : public gsPreconditionerOp<typename MatrixType::Scalar>
279{
280 typedef memory::shared_ptr<MatrixType> MatrixPtr;
281 typedef typename MatrixType::Nested NestedMatrix;
282
283public:
285 typedef typename MatrixType::Scalar T;
286
288 typedef memory::shared_ptr< gsGaussSeidelOp > Ptr;
289
291 typedef memory::unique_ptr< gsGaussSeidelOp > uPtr;
292
295
297 explicit gsGaussSeidelOp(const MatrixType& mat)
298 : m_mat(), m_expr(mat.derived()) {}
299
301 explicit gsGaussSeidelOp(const MatrixPtr& mat)
302 : m_mat(mat), m_expr(m_mat->derived()) { }
303
304 static uPtr make(const MatrixType& mat)
305 { return memory::make_unique( new gsGaussSeidelOp(mat) ); }
306
307 static uPtr make(const MatrixPtr& mat)
308 { return memory::make_unique( new gsGaussSeidelOp(mat) ); }
309
310 void step(const gsMatrix<T> & rhs, gsMatrix<T> & x) const
311 {
312 if ( ordering == gsGaussSeidel::forward )
313 internal::gaussSeidelSweep<T>(m_expr,x,rhs);
314 if ( ordering == gsGaussSeidel::reverse )
315 internal::reverseGaussSeidelSweep<T>(m_expr,x,rhs);
316 if ( ordering == gsGaussSeidel::symmetric )
317 {
318 internal::gaussSeidelSweep<T>(m_expr,x,rhs);
319 internal::reverseGaussSeidelSweep<T>(m_expr,x,rhs);
320 }
321 }
322
323 void stepT(const gsMatrix<T> & rhs, gsMatrix<T> & x) const
324 {
325 if ( ordering == gsGaussSeidel::forward )
326 internal::reverseGaussSeidelSweep<T>(m_expr,x,rhs);
327 if ( ordering == gsGaussSeidel::reverse )
328 internal::gaussSeidelSweep<T>(m_expr,x,rhs);
329 if ( ordering == gsGaussSeidel::symmetric )
330 {
331 internal::gaussSeidelSweep<T>(m_expr,x,rhs);
332 internal::reverseGaussSeidelSweep<T>(m_expr,x,rhs);
333 }
334 }
335
336 index_t rows() const {return m_expr.rows();}
337 index_t cols() const {return m_expr.cols();}
338
340 NestedMatrix matrix() const { return m_expr; }
341
343 MatrixPtr matrixPtr() const {
344 GISMO_ENSURE( m_mat, "A shared pointer is only available if it was provided to gsGaussSeidelOp." );
345 return m_mat;
346 }
347
348 typename gsLinearOperator<T>::Ptr underlyingOp() const { return makeMatrixOp(m_mat); }
349
350private:
351 const MatrixPtr m_mat;
352 NestedMatrix m_expr;
353};
354
357template <class Derived>
358typename gsGaussSeidelOp<Derived>::uPtr makeGaussSeidelOp(const gsEigen::EigenBase<Derived>& mat)
359{ return gsGaussSeidelOp<Derived>::make(mat.derived()); }
360
363template <class Derived>
364typename gsGaussSeidelOp<Derived>::uPtr makeGaussSeidelOp(const memory::shared_ptr<Derived>& mat)
365{ return gsGaussSeidelOp<Derived>::make(mat); }
366
369template <class Derived>
372
375template <class Derived>
378
381template <class Derived>
384
387template <class Derived>
390
395template <typename MatrixType>
396class gsIncompleteLUOp GISMO_FINAL : public gsPreconditionerOp<typename MatrixType::Scalar>
397{
398 typedef memory::shared_ptr<MatrixType> MatrixPtr;
399 typedef typename MatrixType::Nested NestedMatrix;
400
401public:
403 typedef typename MatrixType::Scalar T;
404
406 typedef memory::shared_ptr< gsIncompleteLUOp > Ptr;
407
409 typedef memory::unique_ptr< gsIncompleteLUOp > uPtr;
410
413
415 explicit gsIncompleteLUOp(const MatrixType& mat, index_t fillfactor = 1)
416 : m_mat(), m_expr(mat.derived())
417 {
418 m_ilu.setFillfactor(fillfactor);
419 m_ilu.compute(m_expr);
420 }
421
423 explicit gsIncompleteLUOp(const MatrixPtr& mat, index_t fillfactor = 1)
424 : m_mat(mat), m_expr(m_mat->derived())
425 {
426 m_ilu.setFillfactor(fillfactor);
427 m_ilu.compute(m_expr);
428 }
429
430 static uPtr make(const MatrixType& mat, index_t fillfactor = 1)
431 { return memory::make_unique( new gsIncompleteLUOp(mat,fillfactor) ); }
432
433 static uPtr make(const MatrixPtr& mat, index_t fillfactor = 1)
434 { return memory::make_unique( new gsIncompleteLUOp(mat,fillfactor) ); }
435
436 void step(const gsMatrix<T> & rhs, gsMatrix<T> & x) const
437 {
438 x += m_ilu.solve( rhs - m_expr * x ).eval();
439 }
440
441 // We use our own apply implementation as we can save one multiplication. This is important if the number
442 // of sweeps is 1: Then we can save *all* multiplications.
443 void apply(const gsMatrix<T> & input, gsMatrix<T> & x) const
444 {
445 // For the first sweep, we do not need to multiply with the matrix
446 x = m_ilu.solve( input ).eval();
447
448 for (index_t k = 1; k < Base::m_num_of_sweeps; ++k)
449 x += m_ilu.solve( input - m_expr * x ).eval();
450
451 }
452
453 index_t rows() const {return m_expr.rows();}
454 index_t cols() const {return m_expr.cols();}
455
457 NestedMatrix matrix() const { return m_expr; }
458
460 MatrixPtr matrixPtr() const {
461 GISMO_ENSURE( m_mat, "A shared pointer is only available if it was provided to gsIncompleteLUOp." );
462 return m_mat;
463 }
464
465 typename gsLinearOperator<T>::Ptr underlyingOp() const { return makeMatrixOp(m_mat); }
466
467private:
468 const MatrixPtr m_mat;
469 NestedMatrix m_expr;
470 gsEigen::IncompleteLUT<T> m_ilu;
471 using Base::m_num_of_sweeps;
472};
473
476template <class Derived>
477typename gsIncompleteLUOp<Derived>::uPtr makeIncompleteLUOp(const gsEigen::EigenBase<Derived>& mat)
478{ return gsIncompleteLUOp<Derived>::make(mat.derived()); }
479
482template <class Derived>
483typename gsIncompleteLUOp<Derived>::uPtr makeIncompleteLUOp(const memory::shared_ptr<Derived>& mat)
484{ return gsIncompleteLUOp<Derived>::make(mat); }
485
486
487} // namespace gismo
488
489#ifndef GISMO_BUILD_LIB
490#include GISMO_HPP_HEADER(gsSimplePreconditioners.hpp)
491#endif
Gauss-Seidel preconditioner.
Definition gsSimplePreconditioners.h:279
index_t rows() const
Returns the number of rows of the operator.
Definition gsSimplePreconditioners.h:336
gsGaussSeidelOp< Derived >::uPtr makeGaussSeidelOp(const memory::shared_ptr< Derived > &mat)
Returns a smart pointer to a Jacobi operator referring on mat.
Definition gsSimplePreconditioners.h:364
gsLinearOperator< T >::Ptr underlyingOp() const
Return the underlying operator .
Definition gsSimplePreconditioners.h:348
gsGaussSeidelOp< Derived >::uPtr makeGaussSeidelOp(const gsEigen::EigenBase< Derived > &mat)
Returns a smart pointer to a Gauss-Seidel operator referring on mat.
Definition gsSimplePreconditioners.h:358
const MatrixPtr m_mat
Shared pointer to matrix (if needed)
Definition gsSimplePreconditioners.h:351
gsGaussSeidelOp< Derived, gsGaussSeidel::reverse >::uPtr makeReverseGaussSeidelOp(const memory::shared_ptr< Derived > &mat)
Returns a smart pointer to a reverse Gauss-Seidel operator referring on mat.
Definition gsSimplePreconditioners.h:376
gsGaussSeidelOp(const MatrixPtr &mat)
Constructor with shared pointer to matrix.
Definition gsSimplePreconditioners.h:301
NestedMatrix m_expr
Nested Eigen expression.
Definition gsSimplePreconditioners.h:352
gsGaussSeidelOp< Derived, gsGaussSeidel::symmetric >::uPtr makeSymmetricGaussSeidelOp(const memory::shared_ptr< Derived > &mat)
Returns a smart pointer to a symmetric Gauss-Seidel operator referring on mat.
Definition gsSimplePreconditioners.h:388
MatrixPtr matrixPtr() const
Returns a shared pinter to the matrix.
Definition gsSimplePreconditioners.h:343
memory::shared_ptr< gsGaussSeidelOp > Ptr
Shared pointer for gsGaussSeidelOp.
Definition gsSimplePreconditioners.h:288
NestedMatrix matrix() const
Returns the matrix.
Definition gsSimplePreconditioners.h:340
MatrixType::Scalar T
Scalar type.
Definition gsSimplePreconditioners.h:285
void step(const gsMatrix< T > &rhs, gsMatrix< T > &x) const
Apply the method for given right hand side and current iterate.
Definition gsSimplePreconditioners.h:310
index_t cols() const
Returns the number of columns of the operator.
Definition gsSimplePreconditioners.h:337
gsGaussSeidelOp(const MatrixType &mat)
Constructor with given matrix.
Definition gsSimplePreconditioners.h:297
gsGaussSeidelOp< Derived, gsGaussSeidel::symmetric >::uPtr makeSymmetricGaussSeidelOp(const gsEigen::EigenBase< Derived > &mat)
Returns a smart pointer to a symmetric Gauss-Seidel operator referring on mat.
Definition gsSimplePreconditioners.h:382
memory::unique_ptr< gsGaussSeidelOp > uPtr
Unique pointer for gsGaussSeidelOp.
Definition gsSimplePreconditioners.h:291
gsPreconditionerOp< T > Base
Base class.
Definition gsSimplePreconditioners.h:294
void stepT(const gsMatrix< T > &rhs, gsMatrix< T > &x) const
Apply the transposed variant of the method for given right hand side and current iterate.
Definition gsSimplePreconditioners.h:323
gsGaussSeidelOp< Derived, gsGaussSeidel::reverse >::uPtr makeReverseGaussSeidelOp(const gsEigen::EigenBase< Derived > &mat)
Returns a smart pointer to a reverse Gauss-Seidel operator referring on mat.
Definition gsSimplePreconditioners.h:370
Incomplete LU with thresholding preconditioner.
Definition gsSimplePreconditioners.h:397
index_t rows() const
Returns the number of rows of the operator.
Definition gsSimplePreconditioners.h:453
gsLinearOperator< T >::Ptr underlyingOp() const
Return the underlying operator .
Definition gsSimplePreconditioners.h:465
const MatrixPtr m_mat
Shared pointer to matrix (if needed)
Definition gsSimplePreconditioners.h:468
gsIncompleteLUOp< Derived >::uPtr makeIncompleteLUOp(const gsEigen::EigenBase< Derived > &mat)
Returns a smart pointer to a Gauss-Seidel operator referring on mat.
Definition gsSimplePreconditioners.h:477
memory::unique_ptr< gsIncompleteLUOp > uPtr
Unique pointer for gsIncompleteLUOp.
Definition gsSimplePreconditioners.h:409
memory::shared_ptr< gsIncompleteLUOp > Ptr
Shared pointer for gsIncompleteLUOp.
Definition gsSimplePreconditioners.h:406
NestedMatrix m_expr
Nested Eigen expression.
Definition gsSimplePreconditioners.h:469
MatrixPtr matrixPtr() const
Returns a shared pinter to the matrix.
Definition gsSimplePreconditioners.h:460
gsIncompleteLUOp(const MatrixType &mat, index_t fillfactor=1)
Constructor with given matrix.
Definition gsSimplePreconditioners.h:415
NestedMatrix matrix() const
Returns the matrix.
Definition gsSimplePreconditioners.h:457
MatrixType::Scalar T
Scalar type.
Definition gsSimplePreconditioners.h:403
void step(const gsMatrix< T > &rhs, gsMatrix< T > &x) const
Apply the method for given right hand side and current iterate.
Definition gsSimplePreconditioners.h:436
index_t cols() const
Returns the number of columns of the operator.
Definition gsSimplePreconditioners.h:454
gsIncompleteLUOp(const MatrixPtr &mat, index_t fillfactor=1)
Constructor with shared pointer to matrix.
Definition gsSimplePreconditioners.h:423
gsIncompleteLUOp< Derived >::uPtr makeIncompleteLUOp(const memory::shared_ptr< Derived > &mat)
Returns a smart pointer to a Jacobi operator referring on mat.
Definition gsSimplePreconditioners.h:483
gsPreconditionerOp< T > Base
Base class.
Definition gsSimplePreconditioners.h:412
gsEigen::IncompleteLUT< T > m_ilu
The decomposition itself.
Definition gsSimplePreconditioners.h:470
Jacobi preconditioner.
Definition gsSimplePreconditioners.h:147
index_t rows() const
Returns the number of rows of the operator.
Definition gsSimplePreconditioners.h:204
gsLinearOperator< T >::Ptr underlyingOp() const
Return the underlying operator .
Definition gsSimplePreconditioners.h:237
const MatrixPtr m_mat
Shared pointer to matrix (if needed)
Definition gsSimplePreconditioners.h:240
memory::shared_ptr< gsJacobiOp > Ptr
Shared pointer for gsJacobiOp.
Definition gsSimplePreconditioners.h:156
gsJacobiOp(const MatrixPtr &mat, T tau=1)
Constructor with shared pointer to matrix.
Definition gsSimplePreconditioners.h:169
void setDamping(const T tau)
Set damping parameter.
Definition gsSimplePreconditioners.h:208
NestedMatrix m_expr
Nested Eigen expression.
Definition gsSimplePreconditioners.h:241
MatrixPtr matrixPtr() const
Returns a shared pinter to the matrix.
Definition gsSimplePreconditioners.h:232
void getDamping()
Get damping parameter.
Definition gsSimplePreconditioners.h:211
memory::unique_ptr< gsJacobiOp > uPtr
Unique pointer for gsJacobiOp.
Definition gsSimplePreconditioners.h:159
NestedMatrix matrix() const
Returns the matrix.
Definition gsSimplePreconditioners.h:229
gsJacobiOp< Derived >::uPtr makeJacobiOp(const gsEigen::EigenBase< Derived > &mat, typename Derived::Scalar tau=1)
Returns a smart pointer to a Jacobi operator referring on mat.
Definition gsSimplePreconditioners.h:250
MatrixType::Scalar T
Scalar type.
Definition gsSimplePreconditioners.h:153
void step(const gsMatrix< T > &rhs, gsMatrix< T > &x) const
Apply the method for given right hand side and current iterate.
Definition gsSimplePreconditioners.h:178
virtual void setOptions(const gsOptionList &opt)
Set options based on a gsOptionList object.
Definition gsSimplePreconditioners.h:222
gsJacobiOp(const MatrixType &mat, T tau=1)
Constructor with given matrix.
Definition gsSimplePreconditioners.h:165
index_t cols() const
Returns the number of columns of the operator.
Definition gsSimplePreconditioners.h:205
static gsOptionList defaultOptions()
Get the default options as gsOptionList object.
Definition gsSimplePreconditioners.h:214
gsJacobiOp< Derived >::uPtr makeJacobiOp(const memory::shared_ptr< Derived > &mat, typename Derived::Scalar tau=1)
Returns a smart pointer to a Jacobi operator referring on mat.
Definition gsSimplePreconditioners.h:256
gsPreconditionerOp< T > Base
Base class.
Definition gsSimplePreconditioners.h:162
memory::shared_ptr< gsLinearOperator > Ptr
Shared pointer for gsLinearOperator.
Definition gsLinearOperator.h:33
A matrix with arbitrary coefficient type and fixed or dynamic size.
Definition gsMatrix.h:41
Class which holds a list of parameters/options, and provides easy access to them.
Definition gsOptionList.h:33
Real askReal(const std::string &label, const Real &value=0) const
Reads value for option label from options.
Definition gsOptionList.cpp:139
void addReal(const std::string &label, const std::string &desc, const Real &value)
Adds a option named label, with description desc and value value.
Definition gsOptionList.cpp:211
Simple abstract class for perconditioners.
Definition gsPreconditioner.h:43
virtual void setOptions(const gsOptionList &opt)
Set options based on a gsOptionList object.
Definition gsPreconditioner.h:111
static gsOptionList defaultOptions()
Get the default options as gsOptionList object.
Definition gsPreconditioner.h:103
Richardson preconditioner.
Definition gsSimplePreconditioners.h:34
index_t rows() const
Returns the number of rows of the operator.
Definition gsSimplePreconditioners.h:87
gsRichardsonOp(const MatrixPtr &mat, T tau=1)
Constructor with shared pointer to matrix.
Definition gsSimplePreconditioners.h:56
gsLinearOperator< T >::Ptr underlyingOp() const
Return the underlying operator .
Definition gsSimplePreconditioners.h:120
gsRichardsonOp(const MatrixType &mat, T tau=1)
Constructor with given matrix.
Definition gsSimplePreconditioners.h:52
const MatrixPtr m_mat
Shared pointer to matrix (if needed)
Definition gsSimplePreconditioners.h:123
gsRichardsonOp< Derived >::uPtr makeRichardsonOp(const memory::shared_ptr< Derived > &mat, typename Derived::Scalar tau=1)
Returns a smart pointer to a Richardson operator referring on mat.
Definition gsSimplePreconditioners.h:139
memory::shared_ptr< gsRichardsonOp > Ptr
Shared pointer for gsRichardsonOp.
Definition gsSimplePreconditioners.h:43
gsRichardsonOp< Derived >::uPtr makeRichardsonOp(const gsEigen::EigenBase< Derived > &mat, typename Derived::Scalar tau=1)
Returns a smart pointer to a Richardson operator referring on mat.
Definition gsSimplePreconditioners.h:133
memory::unique_ptr< gsRichardsonOp > uPtr
Unique pointer for gsRichardsonOp.
Definition gsSimplePreconditioners.h:46
void setDamping(const T tau)
Set damping parameter.
Definition gsSimplePreconditioners.h:91
NestedMatrix m_expr
Nested Eigen expression.
Definition gsSimplePreconditioners.h:124
MatrixPtr matrixPtr() const
Returns a shared pinter to the matrix.
Definition gsSimplePreconditioners.h:115
void getDamping()
Get damping parameter.
Definition gsSimplePreconditioners.h:94
NestedMatrix matrix() const
Returns the matrix.
Definition gsSimplePreconditioners.h:112
MatrixType::Scalar T
Scalar type.
Definition gsSimplePreconditioners.h:40
void step(const gsMatrix< T > &rhs, gsMatrix< T > &x) const
Apply the method for given right hand side and current iterate.
Definition gsSimplePreconditioners.h:65
virtual void setOptions(const gsOptionList &opt)
Set options based on a gsOptionList object.
Definition gsSimplePreconditioners.h:105
index_t cols() const
Returns the number of columns of the operator.
Definition gsSimplePreconditioners.h:88
static gsOptionList defaultOptions()
Get the default options as gsOptionList object.
Definition gsSimplePreconditioners.h:97
gsPreconditionerOp< T > Base
Base class.
Definition gsSimplePreconditioners.h:49
ordering
Specify the ordering of gsGaussSeidelOp preconditioner.
Definition gsSimplePreconditioners.h:265
@ symmetric
one total step is = one forward + one backward
Definition gsSimplePreconditioners.h:268
@ reverse
reverse ordering
Definition gsSimplePreconditioners.h:267
@ forward
standard forward ordering
Definition gsSimplePreconditioners.h:266
#define index_t
Definition gsConfig.h:32
#define GISMO_ENSURE(cond, message)
Definition gsDebug.h:102
#define GISMO_ASSERT(cond, message)
Definition gsDebug.h:89
This is the main header file that collects wrappers of Eigen for linear algebra.
Simple abstract class for (discrete) linear preconditioners.
unique_ptr< T > make_unique(T *x)
Definition gsMemory.h:198
The G+Smo namespace, containing all definitions for the library.