4 #ifndef DUNE_SOLVERS_HH
5 #define DUNE_SOLVERS_HH
17 #include <dune/common/timer.hh>
18 #include <dune/common/ftraits.hh>
19 #include <dune/common/static_assert.hh>
93 template<
class X,
class Y>
140 s << std::setw(
normSpacing) <<
"Rate" << std::endl;
144 template <
class DataType>
147 const DataType& norm,
148 const DataType& norm_old)
const
150 const DataType rate = norm/norm_old;
157 template <
class DataType>
160 const DataType& norm)
const
209 template<
class L,
class P>
211 double reduction,
int maxit,
int verbose) :
212 ssp(), _op(op), _prec(prec), _sp(ssp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
214 dune_static_assert(static_cast<int>(L::category) == static_cast<int>(P::category),
215 "L and P have to have the same category!");
217 "L has to be sequential!");
240 template<
class L,
class S,
class P>
242 double reduction,
int maxit,
int verbose) :
243 _op(op), _prec(prec), _sp(sp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
245 dune_static_assert( static_cast<int>(L::category) == static_cast<int>(P::category),
246 "L and P must have the same category!");
247 dune_static_assert( static_cast<int>(L::category) == static_cast<int>(S::category),
248 "L and S must have the same category!");
268 double def0 = _sp.norm(b);
273 std::cout <<
"=== LoopSolver" << std::endl;
285 int i=1;
double def=def0;
286 for ( ; i<=_maxit; i++ )
292 double defnew=_sp.norm(b);
297 if (def<def0*_reduction || def<1E-30)
320 std::cout <<
"=== rate=" << res.
conv_rate
323 <<
", IT=" << i << std::endl;
330 std::swap(_reduction,reduction);
331 (*this).apply(x,b,res);
332 std::swap(_reduction,reduction);
363 template<
class L,
class P>
365 double reduction,
int maxit,
int verbose) :
366 ssp(), _op(op), _prec(prec), _sp(ssp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
368 dune_static_assert(static_cast<int>(L::category) == static_cast<int>(P::category),
369 "L and P have to have the same category!");
371 "L has to be sequential!");
378 template<
class L,
class S,
class P>
380 double reduction,
int maxit,
int verbose) :
381 _op(op), _prec(prec), _sp(sp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
383 dune_static_assert(static_cast<int>(L::category) == static_cast<int>(P::category),
384 "L and P have to have the same category!");
385 dune_static_assert(static_cast<int>(L::category) == static_cast<int>(S::category),
386 "L and S have to have the same category!");
404 double def0 = _sp.norm(b);
408 std::cout <<
"=== GradientSolver" << std::endl;
416 int i=1;
double def=def0;
418 for ( ; i<=_maxit; i++ )
423 lambda = _sp.dot(p,b)/_sp.dot(q,p);
427 double defnew=_sp.norm(b);
432 if (def<def0*_reduction || def<1E-30)
448 std::cout <<
"=== rate=" << res.
conv_rate
451 <<
", IT=" << i << std::endl;
461 std::swap(_reduction,reduction);
462 (*this).apply(x,b,res);
463 std::swap(_reduction,reduction);
494 template<
class L,
class P>
495 CGSolver (L& op, P& prec,
double reduction,
int maxit,
int verbose) :
496 ssp(), _op(op), _prec(prec), _sp(ssp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
498 dune_static_assert( static_cast<int>(L::category) == static_cast<int>(P::category),
499 "L and P must have the same category!");
501 "L must be sequential!");
508 template<
class L,
class S,
class P>
509 CGSolver (L& op, S& sp, P& prec,
double reduction,
int maxit,
int verbose) :
510 _op(op), _prec(prec), _sp(sp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
512 dune_static_assert( static_cast<int>(L::category) == static_cast<int>(P::category),
513 "L and P must have the same category!");
514 dune_static_assert( static_cast<int>(L::category) == static_cast<int>(S::category),
515 "L and S must have the same category!");
533 double def0 = _sp.norm(b);
542 std::cout <<
"=== rate=" << res.
conv_rate
544 <<
", IT=0" << std::endl;
550 std::cout <<
"=== CGSolver" << std::endl;
564 rholast = _sp.dot(p,b);
568 for ( ; i<=_maxit; i++ )
572 alpha = _sp.dot(p,q);
573 lambda = rholast/alpha;
578 double defnew=_sp.norm(b);
584 if (def<def0*_reduction || def<1E-30)
611 std::cout <<
"=== rate=" << res.
conv_rate
614 <<
", IT=" << i << std::endl;
623 virtual void apply (X& x, X& b,
double reduction,
626 std::swap(_reduction,reduction);
627 (*this).apply(x,b,res);
628 std::swap(_reduction,reduction);
659 template<
class L,
class P>
661 double reduction,
int maxit,
int verbose) :
662 ssp(), _op(op), _prec(prec), _sp(ssp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
664 dune_static_assert(static_cast<int>(L::category) == static_cast<int>(P::category),
"L and P must be of the same category!");
672 template<
class L,
class S,
class P>
674 double reduction,
int maxit,
int verbose) :
675 _op(op), _prec(prec), _sp(sp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
677 dune_static_assert( static_cast<int>(L::category) == static_cast<int>(P::category),
678 "L and P must have the same category!");
679 dune_static_assert( static_cast<int>(L::category) == static_cast<int>(S::category),
680 "L and S must have the same category!");
690 const double EPSILON=1e-80;
693 field_type rho, rho_new, alpha, beta, h, omega;
718 norm = norm_old = norm_0 = _sp.norm(r);
729 std::cout <<
"=== BiCGSTABSolver" << std::endl;
739 if ( norm < (_reduction * norm_0) || norm<1E-30)
754 for (it = 0.5; it < _maxit; it+=.5)
761 rho_new = _sp.dot(rt,r);
764 if (std::abs(rho) <= EPSILON)
765 DUNE_THROW(
ISTLError,
"breakdown in BiCGSTAB - rho "
766 << rho <<
" <= EPSILON " << EPSILON
767 <<
" after " << it <<
" iterations");
768 if (std::abs(omega) <= EPSILON)
769 DUNE_THROW(
ISTLError,
"breakdown in BiCGSTAB - omega "
770 << omega <<
" <= EPSILON " << EPSILON
771 <<
" after " << it <<
" iterations");
778 beta = ( rho_new / rho ) * ( alpha / omega );
794 if ( std::abs(h) < EPSILON )
817 if ( norm < (_reduction * norm_0) )
834 omega = _sp.dot(t,r)/_sp.dot(t,t);
856 if ( norm < (_reduction * norm_0) || norm<1E-30)
869 res.
iterations =
static_cast<int>(std::ceil(it));
874 std::cout <<
"=== rate=" << res.
conv_rate
877 <<
", IT=" << it << std::endl;
887 std::swap(_reduction,reduction);
888 (*this).apply(x,b,res);
889 std::swap(_reduction,reduction);
918 typedef typename FieldTraits<field_type>::real_type
real_type;
925 template<
class L,
class P>
926 MINRESSolver (L& op, P& prec,
double reduction,
int maxit,
int verbose) :
927 ssp(), _op(op), _prec(prec), _sp(ssp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
929 dune_static_assert( static_cast<int>(L::category) == static_cast<int>(P::category),
930 "L and P must have the same category!");
932 "L must be sequential!");
939 template<
class L,
class S,
class P>
940 MINRESSolver (L& op, S& sp, P& prec,
double reduction,
int maxit,
int verbose) :
941 _op(op), _prec(prec), _sp(sp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
943 dune_static_assert( static_cast<int>(L::category) == static_cast<int>(P::category),
944 "L and P must have the same category!");
945 dune_static_assert( static_cast<int>(L::category) == static_cast<int>(S::category),
946 "L and S must have the same category!");
971 std::cout <<
"=== rate=" << res.
conv_rate <<
", T=" << res.
elapsed <<
", TIT=" << res.
elapsed <<
", IT=0" << std::endl;
977 std::cout <<
"=== MINRESSolver" << std::endl;
1003 beta = std::sqrt(std::abs(_sp.dot(z,b)));
1009 q[0].resize(b.size());
1010 q[1].resize(b.size());
1011 q[2].resize(b.size());
1017 p[0].resize(b.size());
1018 p[1].resize(b.size());
1019 p[2].resize(b.size());
1029 for ( ; i<=_maxit; i++)
1039 q[i2].axpy(-beta, q[i0]);
1040 alpha = _sp.dot(q[i2],z);
1041 q[i2].axpy(-alpha, q[i1]);
1044 _prec.
apply(z,q[i2]);
1046 beta = std::sqrt(std::abs(_sp.dot(q[i2],z)));
1061 T[2] = c[(i+1)%2]*alpha - s[(i+1)%2]*T[1];
1062 T[1] = c[(i+1)%2]*T[1] + s[(i+1)%2]*alpha;
1069 c[i%2] = 1.0/std::sqrt(T[2]*T[2] + beta*beta);
1070 s[i%2] = beta*c[i%2];
1074 T[2] = c[i%2]*T[2] + s[i%2]*beta;
1077 xi[i%2] = -s[i%2]*xi[(i+1)%2];
1078 xi[(i+1)%2] *= c[i%2];
1082 p[i2].axpy(-T[1],p[i1]);
1083 p[i2].axpy(-T[0],p[i0]);
1087 x.axpy(beta0*xi[(i+1)%2], p[i2]);
1098 real_type defnew = std::abs(beta0*xi[i%2]);
1104 if (def<def0*_reduction || def<1E-30 || i==_maxit)
1118 res.
elapsed = watch.elapsed();
1122 std::cout <<
"=== rate=" << res.
conv_rate
1125 <<
", IT=" << i << std::endl;
1137 std::swap(_reduction,reduction);
1138 (*this).apply(x,b,res);
1139 std::swap(_reduction,reduction);
1163 template<
class X,
class Y=X,
class F = Y>
1183 template<
class L,
class P>
1184 RestartedGMResSolver (L& op, P& prec,
double reduction,
int restart,
int maxit,
int verbose,
bool recalc_defect =
false) :
1186 ssp(), _sp(ssp), _restart(restart),
1187 _reduction(reduction), _maxit(maxit), _verbose(verbose),
1188 _recalc_defect(recalc_defect)
1190 dune_static_assert(static_cast<int>(P::category) == static_cast<int>(L::category),
1191 "P and L must be the same category!");
1193 "L must be sequential!");
1203 template<
class L,
class S,
class P>
1204 RestartedGMResSolver (L& op, S& sp, P& prec,
double reduction,
int restart,
int maxit,
int verbose,
bool recalc_defect =
false) :
1206 _sp(sp), _restart(restart),
1207 _reduction(reduction), _maxit(maxit), _verbose(verbose),
1208 _recalc_defect(recalc_defect)
1210 dune_static_assert(static_cast<int>(P::category) == static_cast<int>(L::category),
1211 "P and L must have the same category!");
1212 dune_static_assert(static_cast<int>(P::category) == static_cast<int>(S::category),
1213 "P and S must have the same category!");
1219 apply(x,b,_reduction,res);
1235 std::vector<field_type> s(m+1), cs(m), sn(m);
1238 std::vector< std::vector<field_type> > H(m+1,s);
1239 std::vector<F> v(m+1,b);
1250 w = 0.0; _M.
apply(w,b);
1251 norm_0 = _sp.norm(w);
1255 v[0] = 0.0; _M.
apply(v[0],w);
1256 beta = _sp.norm(v[0]);
1261 w = 0.0; _M.
apply(w,b);
1264 norm_0 = _sp.norm(b);
1265 v[0] = 0.0; _M.
apply(v[0],b);
1266 beta = _sp.norm(v[0]);
1272 norm = norm_old = _sp.norm(v[0]);
1277 std::cout <<
"=== RestartedGMResSolver" << std::endl;
1287 if (norm <= reduction * norm_0) {
1295 while (j <= _maxit && res.
converged !=
true) {
1296 v[0] *= (1.0 / beta);
1297 for (i=1; i<=m; i++) s[i] = 0.0;
1300 for (i = 0; i < m && j <= _maxit && res.
converged !=
true; i++, j++) {
1303 _A_.
apply(v[i], v[i+1]);
1304 _M.
apply(w, v[i+1]);
1305 for (k = 0; k <= i; k++) {
1306 H[k][i] = _sp.dot(w, v[k]);
1308 w.axpy(-H[k][i], v[k]);
1310 H[i+1][i] = _sp.norm(w);
1311 if (H[i+1][i] == 0.0)
1312 DUNE_THROW(
ISTLError,
"breakdown in GMRes - |w| "
1313 << w <<
" == 0.0 after " << j <<
" iterations");
1315 v[i+1] = w; v[i+1] *= (1.0 / H[i+1][i]);
1317 for (k = 0; k < i; k++)
1318 applyPlaneRotation(H[k][i], H[k+1][i], cs[k], sn[k]);
1320 generatePlaneRotation(H[i][i], H[i+1][i], cs[i], sn[i]);
1321 applyPlaneRotation(H[i][i], H[i+1][i], cs[i], sn[i]);
1322 applyPlaneRotation(s[i], s[i+1], cs[i], sn[i]);
1324 norm = std::abs(s[i+1]);
1333 if (norm < reduction * norm_0) {
1341 update(x, i - 1, H, s, v);
1347 beta = _sp.norm(v[0]);
1354 update(w, i - 1, H, s, v);
1363 v[0] = 0.0; _M.
apply(v[0],b);
1364 beta = _sp.norm(v[0]);
1377 if (norm < reduction * norm_0) {
1382 if (res.
converged !=
true && _verbose > 0)
1383 std::cout <<
"=== GMRes::restart\n";
1391 res.
elapsed = watch.elapsed();
1402 std::cout <<
"=== rate=" << res.
conv_rate
1411 std::vector< std::vector<field_type> > & h,
1412 std::vector<field_type> & s, std::vector<F> v)
1414 std::vector<field_type> y(s);
1417 for (
int i = k; i >= 0; i--) {
1419 for (
int j = i - 1; j >= 0; j--)
1420 y[j] -= h[j][i] * y[i];
1423 for (
int j = 0; j <= k; j++)
1434 }
else if (std::abs(dy) > std::abs(dx)) {
1436 sn = 1.0 / std::sqrt( 1.0 + temp*temp );
1440 cs = 1.0 / std::sqrt( 1.0 + temp*temp );
1450 dy = -sn * dx + cs * dy;
1454 LinearOperator<X,X>& _A_;
1455 Preconditioner<X,X>& _M;
1456 SeqScalarProduct<X> ssp;
1457 ScalarProduct<X>& _sp;
1462 bool _recalc_defect;