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);
654 typedef typename FieldTraits<field_type>::real_type
real_type;
661 template<
class L,
class P>
663 double reduction,
int maxit,
int verbose) :
664 ssp(), _op(op), _prec(prec), _sp(ssp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
666 dune_static_assert(static_cast<int>(L::category) == static_cast<int>(P::category),
"L and P must be of the same category!");
674 template<
class L,
class S,
class P>
676 double reduction,
int maxit,
int verbose) :
677 _op(op), _prec(prec), _sp(sp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
679 dune_static_assert( static_cast<int>(L::category) == static_cast<int>(P::category),
680 "L and P must have the same category!");
681 dune_static_assert( static_cast<int>(L::category) == static_cast<int>(S::category),
682 "L and S must have the same category!");
692 const double EPSILON=1e-80;
694 field_type rho, rho_new, alpha, beta, h, omega;
719 norm = norm_old = norm_0 = _sp.norm(r);
730 std::cout <<
"=== BiCGSTABSolver" << std::endl;
740 if ( norm < (_reduction * norm_0) || norm<1E-30)
755 for (it = 0.5; it < _maxit; it+=.5)
762 rho_new = _sp.dot(rt,r);
765 if (std::abs(rho) <= EPSILON)
766 DUNE_THROW(
ISTLError,
"breakdown in BiCGSTAB - rho "
767 << rho <<
" <= EPSILON " << EPSILON
768 <<
" after " << it <<
" iterations");
769 if (std::abs(omega) <= EPSILON)
770 DUNE_THROW(
ISTLError,
"breakdown in BiCGSTAB - omega "
771 << omega <<
" <= EPSILON " << EPSILON
772 <<
" after " << it <<
" iterations");
779 beta = ( rho_new / rho ) * ( alpha / omega );
795 if ( std::abs(h) < EPSILON )
818 if ( norm < (_reduction * norm_0) )
835 omega = _sp.dot(t,r)/_sp.dot(t,t);
857 if ( norm < (_reduction * norm_0) || norm<1E-30)
870 res.
iterations =
static_cast<int>(std::ceil(it));
875 std::cout <<
"=== rate=" << res.
conv_rate
878 <<
", IT=" << it << std::endl;
888 std::swap(_reduction,reduction);
889 (*this).apply(x,b,res);
890 std::swap(_reduction,reduction);
919 typedef typename FieldTraits<field_type>::real_type
real_type;
926 template<
class L,
class P>
927 MINRESSolver (L& op, P& prec,
double reduction,
int maxit,
int verbose) :
928 ssp(), _op(op), _prec(prec), _sp(ssp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
930 dune_static_assert( static_cast<int>(L::category) == static_cast<int>(P::category),
931 "L and P must have the same category!");
933 "L must be sequential!");
940 template<
class L,
class S,
class P>
941 MINRESSolver (L& op, S& sp, P& prec,
double reduction,
int maxit,
int verbose) :
942 _op(op), _prec(prec), _sp(sp), _reduction(reduction), _maxit(maxit), _verbose(verbose)
944 dune_static_assert( static_cast<int>(L::category) == static_cast<int>(P::category),
945 "L and P must have the same category!");
946 dune_static_assert( static_cast<int>(L::category) == static_cast<int>(S::category),
947 "L and S must have the same category!");
972 std::cout <<
"=== rate=" << res.
conv_rate <<
", T=" << res.
elapsed <<
", TIT=" << res.
elapsed <<
", IT=0" << std::endl;
978 std::cout <<
"=== MINRESSolver" << std::endl;
1004 beta = std::sqrt(std::abs(_sp.dot(z,b)));
1010 q[0].resize(b.size());
1011 q[1].resize(b.size());
1012 q[2].resize(b.size());
1018 p[0].resize(b.size());
1019 p[1].resize(b.size());
1020 p[2].resize(b.size());
1030 for ( ; i<=_maxit; i++)
1040 q[i2].axpy(-beta, q[i0]);
1041 alpha = _sp.dot(q[i2],z);
1042 q[i2].axpy(-alpha, q[i1]);
1045 _prec.
apply(z,q[i2]);
1047 beta = std::sqrt(std::abs(_sp.dot(q[i2],z)));
1062 T[2] = c[(i+1)%2]*alpha - s[(i+1)%2]*T[1];
1063 T[1] = c[(i+1)%2]*T[1] + s[(i+1)%2]*alpha;
1070 c[i%2] = 1.0/std::sqrt(T[2]*T[2] + beta*beta);
1071 s[i%2] = beta*c[i%2];
1075 T[2] = c[i%2]*T[2] + s[i%2]*beta;
1078 xi[i%2] = -s[i%2]*xi[(i+1)%2];
1079 xi[(i+1)%2] *= c[i%2];
1083 p[i2].axpy(-T[1],p[i1]);
1084 p[i2].axpy(-T[0],p[i0]);
1088 x.axpy(beta0*xi[(i+1)%2], p[i2]);
1099 real_type defnew = std::abs(beta0*xi[i%2]);
1105 if (def<def0*_reduction || def<1E-30 || i==_maxit)
1119 res.
elapsed = watch.elapsed();
1123 std::cout <<
"=== rate=" << res.
conv_rate
1126 <<
", IT=" << i << std::endl;
1138 std::swap(_reduction,reduction);
1139 (*this).apply(x,b,res);
1140 std::swap(_reduction,reduction);
1164 template<
class X,
class Y=X,
class F = Y>
1175 typedef typename FieldTraits<field_type>::real_type
real_type;
1186 template<
class L,
class P>
1187 RestartedGMResSolver (L& op, P& prec,
double reduction,
int restart,
int maxit,
int verbose,
bool recalc_defect =
false) :
1189 ssp(), _sp(ssp), _restart(restart),
1190 _reduction(reduction), _maxit(maxit), _verbose(verbose),
1191 _recalc_defect(recalc_defect)
1193 dune_static_assert(static_cast<int>(P::category) == static_cast<int>(L::category),
1194 "P and L must be the same category!");
1196 "L must be sequential!");
1206 template<
class L,
class S,
class P>
1207 RestartedGMResSolver (L& op, S& sp, P& prec,
double reduction,
int restart,
int maxit,
int verbose,
bool recalc_defect =
false) :
1209 _sp(sp), _restart(restart),
1210 _reduction(reduction), _maxit(maxit), _verbose(verbose),
1211 _recalc_defect(recalc_defect)
1213 dune_static_assert(static_cast<int>(P::category) == static_cast<int>(L::category),
1214 "P and L must have the same category!");
1215 dune_static_assert(static_cast<int>(P::category) == static_cast<int>(S::category),
1216 "P and S must have the same category!");
1222 apply(x,b,_reduction,res);
1238 std::vector<field_type> s(m+1), cs(m), sn(m);
1241 std::vector< std::vector<field_type> > H(m+1,s);
1242 std::vector<F> v(m+1,b);
1253 w = 0.0; _M.
apply(w,b);
1254 norm_0 = _sp.norm(w);
1258 v[0] = 0.0; _M.
apply(v[0],w);
1259 beta = _sp.norm(v[0]);
1265 norm_0 = _sp.norm(b);
1266 v[0] = 0.0; _M.
apply(v[0],b);
1267 beta = _sp.norm(v[0]);
1273 norm = norm_old = _sp.norm(v[0]);
1278 std::cout <<
"=== RestartedGMResSolver" << std::endl;
1288 if (norm <= reduction * norm_0) {
1296 while (j <= _maxit && res.
converged !=
true) {
1297 v[0] *= (1.0 / beta);
1298 for (i=1; i<=m; i++) s[i] = 0.0;
1301 for (i = 0; i < m && j <= _maxit && res.
converged !=
true; i++, j++) {
1304 _A_.
apply(v[i], v[i+1]);
1305 _M.
apply(w, v[i+1]);
1306 for (k = 0; k <= i; k++) {
1307 H[k][i] = _sp.dot(w, v[k]);
1309 w.axpy(-H[k][i], v[k]);
1311 H[i+1][i] = _sp.norm(w);
1312 if (H[i+1][i] == 0.0)
1313 DUNE_THROW(
ISTLError,
"breakdown in GMRes - |w| "
1314 << w <<
" == 0.0 after " << j <<
" iterations");
1316 v[i+1] = w; v[i+1] *= (1.0 / H[i+1][i]);
1318 for (k = 0; k < i; k++)
1319 applyPlaneRotation(H[k][i], H[k+1][i], cs[k], sn[k]);
1321 generatePlaneRotation(H[i][i], H[i+1][i], cs[i], sn[i]);
1322 applyPlaneRotation(H[i][i], H[i+1][i], cs[i], sn[i]);
1323 applyPlaneRotation(s[i], s[i+1], cs[i], sn[i]);
1325 norm = std::abs(s[i+1]);
1334 if (norm < reduction * norm_0) {
1342 update(x, i - 1, H, s, v);
1348 beta = _sp.norm(v[0]);
1355 update(w, i - 1, H, s, v);
1364 v[0] = 0.0; _M.
apply(v[0],b);
1365 beta = _sp.norm(v[0]);
1378 if (norm < reduction * norm_0) {
1383 if (res.
converged !=
true && _verbose > 0)
1384 std::cout <<
"=== GMRes::restart\n";
1392 res.
elapsed = watch.elapsed();
1403 std::cout <<
"=== rate=" << res.
conv_rate
1412 std::vector< std::vector<field_type> > & h,
1413 std::vector<field_type> & s, std::vector<F> v)
1415 std::vector<field_type> y(s);
1418 for (
int i = k; i >= 0; i--) {
1420 for (
int j = i - 1; j >= 0; j--)
1421 y[j] -= h[j][i] * y[i];
1424 for (
int j = 0; j <= k; j++)
1435 }
else if (std::abs(dy) > std::abs(dx)) {
1437 sn = 1.0 / std::sqrt( 1.0 + temp*temp );
1441 cs = 1.0 / std::sqrt( 1.0 + temp*temp );
1451 dy = -sn * dx + cs * dy;
1455 LinearOperator<X,X>& _A_;
1456 Preconditioner<X,X>& _M;
1457 SeqScalarProduct<X> ssp;
1458 ScalarProduct<X>& _sp;
1463 bool _recalc_defect;