39 #ifndef GETFEM_MODEL_SOLVERS_H__
40 #define GETFEM_MODEL_SOLVERS_H__
53 template <
typename T>
struct sort_abs_val_
54 {
bool operator()(T x, T y) {
return (gmm::abs(x) < gmm::abs(y)); } };
56 template <
typename MAT>
void print_eigval(
const MAT &M) {
58 typedef typename gmm::linalg_traits<MAT>::value_type T;
60 gmm::dense_matrix<T> MM(n, n), Q(n, n);
61 std::vector<T> eigval(n);
64 gmm::implicit_qr_algorithm(MM, eigval, Q);
65 std::sort(eigval.begin(), eigval.end(), sort_abs_val_<T>());
66 cout <<
"eival = " << eigval << endl;
80 template <
typename MAT,
typename VECT>
81 struct abstract_linear_solver {
84 virtual void operator ()(
const MAT &, VECT &,
const VECT &,
86 virtual ~abstract_linear_solver() {}
89 template <
typename MAT,
typename VECT>
90 struct linear_solver_cg_preconditioned_ildlt
91 :
public abstract_linear_solver<MAT, VECT> {
92 void operator ()(
const MAT &M, VECT &x,
const VECT &b,
95 gmm::cg(M, x, b, P, iter);
96 if (!iter.converged()) GMM_WARNING2(
"cg did not converge!");
100 template <
typename MAT,
typename VECT>
101 struct linear_solver_gmres_preconditioned_ilu
102 :
public abstract_linear_solver<MAT, VECT> {
103 void operator ()(
const MAT &M, VECT &x,
const VECT &b,
107 if (!iter.converged()) GMM_WARNING2(
"gmres did not converge!");
111 template <
typename MAT,
typename VECT>
112 struct linear_solver_gmres_unpreconditioned
113 :
public abstract_linear_solver<MAT, VECT> {
114 void operator ()(
const MAT &M, VECT &x,
const VECT &b,
116 gmm::identity_matrix P;
118 if (!iter.converged()) GMM_WARNING2(
"gmres did not converge!");
122 template <
typename MAT,
typename VECT>
123 struct linear_solver_gmres_preconditioned_ilut
124 :
public abstract_linear_solver<MAT, VECT> {
125 void operator ()(
const MAT &M, VECT &x,
const VECT &b,
129 if (!iter.converged()) GMM_WARNING2(
"gmres did not converge!");
133 template <
typename MAT,
typename VECT>
134 struct linear_solver_gmres_preconditioned_ilutp
135 :
public abstract_linear_solver<MAT, VECT> {
136 void operator ()(
const MAT &M, VECT &x,
const VECT &b,
140 if (!iter.converged()) GMM_WARNING2(
"gmres did not converge!");
144 #if defined(GMM_USES_SUPERLU)
145 template <
typename MAT,
typename VECT>
146 struct linear_solver_superlu
147 :
public abstract_linear_solver<MAT, VECT> {
148 void operator ()(
const MAT &M, VECT &x,
const VECT &b,
154 int info = gmm::SuperLU_solve(M, x, b, rcond);
155 iter.enforce_converged(info == 0);
156 if (iter.get_noisy()) cout <<
"condition number: " << 1.0/rcond<< endl;
161 template <
typename MAT,
typename VECT>
162 struct linear_solver_dense_lu :
public abstract_linear_solver<MAT, VECT> {
163 void operator ()(
const MAT &M, VECT &x,
const VECT &b,
165 typedef typename gmm::linalg_traits<MAT>::value_type T;
166 gmm::dense_matrix<T> MM(gmm::mat_nrows(M),gmm::mat_ncols(M));
169 iter.enforce_converged(
true);
173 #if defined(GMM_USES_MUMPS)
174 template <
typename MAT,
typename VECT>
175 struct linear_solver_mumps :
public abstract_linear_solver<MAT, VECT> {
176 void operator ()(
const MAT &M, VECT &x,
const VECT &b,
178 bool ok = gmm::MUMPS_solve(M, x, b,
false);
179 iter.enforce_converged(ok);
182 template <
typename MAT,
typename VECT>
183 struct linear_solver_mumps_sym :
public abstract_linear_solver<MAT, VECT> {
184 void operator ()(
const MAT &M, VECT &x,
const VECT &b,
186 bool ok = gmm::MUMPS_solve(M, x, b,
true);
187 iter.enforce_converged(ok);
192 #if GETFEM_PARA_LEVEL > 1 && GETFEM_PARA_SOLVER == MUMPS_PARA_SOLVER
193 template <
typename MAT,
typename VECT>
194 struct linear_solver_distributed_mumps
195 :
public abstract_linear_solver<MAT, VECT> {
196 void operator ()(
const MAT &M, VECT &x,
const VECT &b,
198 double tt_ref=MPI_Wtime();
199 bool ok = MUMPS_distributed_matrix_solve(M, x, b,
false);
200 iter.enforce_converged(ok);
201 if (MPI_IS_MASTER()) cout<<
"UNSYMMETRIC MUMPS time "<< MPI_Wtime() - tt_ref<<endl;
205 template <
typename MAT,
typename VECT>
206 struct linear_solver_distributed_mumps_sym
207 :
public abstract_linear_solver<MAT, VECT> {
208 void operator ()(
const MAT &M, VECT &x,
const VECT &b,
210 double tt_ref=MPI_Wtime();
211 bool ok = MUMPS_distributed_matrix_solve(M, x, b,
true);
212 iter.enforce_converged(ok);
213 if (MPI_IS_MASTER()) cout<<
"SYMMETRIC MUMPS time "<< MPI_Wtime() - tt_ref<<endl;
223 struct abstract_newton_line_search {
224 double conv_alpha, conv_r;
225 size_t it, itmax, glob_it;
227 virtual void init_search(
double r,
size_t git,
double R0 = 0.0) = 0;
228 virtual double next_try(
void) = 0;
229 virtual bool is_converged(
double,
double R1 = 0.0) = 0;
230 virtual double converged_value(
void) {
234 virtual double converged_residual(
void) {
return conv_r; };
235 virtual ~abstract_newton_line_search() { }
239 struct newton_search_with_step_control :
public abstract_newton_line_search {
241 virtual void init_search(
double ,
size_t ,
double = 0.0)
242 { GMM_ASSERT1(
false,
"Not to be used"); }
244 virtual double next_try(
void)
245 { GMM_ASSERT1(
false,
"Not to be used"); }
247 virtual bool is_converged(
double ,
double = 0.0)
248 { GMM_ASSERT1(
false,
"Not to be used"); }
250 newton_search_with_step_control() {}
254 struct quadratic_newton_line_search :
public abstract_newton_line_search {
257 double alpha, first_res;
258 virtual void init_search(
double r,
size_t git,
double R0 = 0.0) {
259 GMM_ASSERT1(R0 != 0.0,
"You have to specify R0");
261 conv_alpha =
alpha = double(1); conv_r = first_res = r; it = 0;
264 virtual double next_try(
void) {
266 if (it == 1)
return double(1);
267 GMM_ASSERT1(R1_ != 0.0,
"You have to specify R1");
268 double a = R0_ / R1_;
269 return (a < 0) ? (a*0.5 + sqrt(a*a*0.25-a)) : a*0.5;
271 virtual bool is_converged(
double r,
double R1 = 0.0) {
273 R1_ = R1;
return ((gmm::abs(R1_) < gmm::abs(R0_*0.5)) || it >= itmax);
275 quadratic_newton_line_search(
size_t imax =
size_t(-1)) { itmax = imax; }
279 struct simplest_newton_line_search :
public abstract_newton_line_search {
280 double alpha, alpha_mult, first_res, alpha_max_ratio, alpha_min,
282 virtual void init_search(
double r,
size_t git,
double = 0.0) {
284 conv_alpha =
alpha = double(1); conv_r = first_res = r; it = 0;
286 virtual double next_try(
void)
287 { conv_alpha =
alpha;
alpha *= alpha_mult; ++it;
return conv_alpha; }
288 virtual bool is_converged(
double r,
double = 0.0) {
290 return ((it <= 1 && r < first_res)
291 || (r <= first_res * alpha_max_ratio && r <= alpha_threshold_res)
292 || (conv_alpha <= alpha_min && r < first_res * 1e5)
295 simplest_newton_line_search
296 (
size_t imax =
size_t(-1),
double a_max_ratio = 6.0/5.0,
297 double a_min = 1.0/1000.0,
double a_mult = 3.0/5.0,
298 double a_threshold_res = 1e50)
299 : alpha_mult(a_mult), alpha_max_ratio(a_max_ratio), alpha_min(a_min),
300 alpha_threshold_res(a_threshold_res)
304 struct default_newton_line_search :
public abstract_newton_line_search {
322 double alpha, alpha_old, alpha_mult, first_res, alpha_max_ratio;
323 double alpha_min_ratio, alpha_min;
325 bool max_ratio_reached;
326 double alpha_max_ratio_reached, r_max_ratio_reached;
329 virtual void init_search(
double r,
size_t git,
double = 0.0);
330 virtual double next_try(
void);
331 virtual bool is_converged(
double r,
double = 0.0);
332 default_newton_line_search(
void) { count_pat = 0; }
336 struct basic_newton_line_search :
public abstract_newton_line_search {
337 double alpha, alpha_mult, first_res, alpha_max_ratio;
338 double alpha_min, prev_res, alpha_max_augment;
339 virtual void init_search(
double r,
size_t git,
double = 0.0) {
341 conv_alpha =
alpha = double(1);
342 prev_res = conv_r = first_res = r; it = 0;
344 virtual double next_try(
void)
345 { conv_alpha =
alpha;
alpha *= alpha_mult; ++it;
return conv_alpha; }
346 virtual bool is_converged(
double r,
double = 0.0) {
347 if (glob_it == 0 || (r < first_res /
double(2))
348 || (conv_alpha <= alpha_min && r < first_res * alpha_max_augment)
350 { conv_r = r;
return true; }
351 if (it > 1 && r > prev_res && prev_res < alpha_max_ratio * first_res)
353 prev_res = conv_r = r;
356 basic_newton_line_search
357 (
size_t imax =
size_t(-1),
358 double a_max_ratio = 5.0/3.0,
359 double a_min = 1.0/1000.0,
double a_mult = 3.0/5.0,
double a_augm = 2.0)
360 : alpha_mult(a_mult), alpha_max_ratio(a_max_ratio),
361 alpha_min(a_min), alpha_max_augment(a_augm) { itmax = imax; }
365 struct systematic_newton_line_search :
public abstract_newton_line_search {
366 double alpha, alpha_mult, first_res;
367 double alpha_min, prev_res;
369 virtual void init_search(
double r,
size_t git,
double = 0.0) {
371 conv_alpha =
alpha = double(1);
372 prev_res = conv_r = first_res = r; it = 0; first =
true;
374 virtual double next_try(
void)
375 {
double a =
alpha;
alpha *= alpha_mult; ++it;
return a; }
376 virtual bool is_converged(
double r,
double = 0.0) {
378 if (r < conv_r || first)
379 { conv_r = r; conv_alpha =
alpha / alpha_mult; first =
false; }
380 if ((alpha <= alpha_min*alpha_mult) || it >= itmax)
return true;
383 systematic_newton_line_search
384 (
size_t imax =
size_t(-1),
385 double a_min = 1.0/10000.0,
double a_mult = 3.0/5.0)
386 : alpha_mult(a_mult), alpha_min(a_min) { itmax = imax; }
414 template <
typename PB>
417 typedef typename gmm::linalg_traits<typename PB::VECTOR>::value_type T;
418 typedef typename gmm::number_traits<T>::magnitude_type R;
420 iter_linsolv0.reduce_noisy();
421 iter_linsolv0.set_resmax(iter.get_resmax()/20.0);
422 iter_linsolv0.set_maxiter(10000);
424 pb.compute_residual();
425 R approx_eln = pb.approx_external_load_norm();
426 if (approx_eln == R(0)) approx_eln = R(1);
428 typename PB::VECTOR b0(gmm::vect_size(pb.residual()));
430 typename PB::VECTOR Xi(gmm::vect_size(pb.residual()));
433 typename PB::VECTOR dr(gmm::vect_size(pb.residual()));
435 R alpha0(0),
alpha(1), res0(gmm::vect_norm1(b0)), minres(res0);
442 scalar_type crit = pb.residual_norm() / approx_eln;
443 if (iter.finished(crit))
return;
448 if (!iter.converged(crit)) {
452 while (is_singular) {
454 pb.add_to_residual(b0, alpha-R(1));
456 if (iter.get_noisy() > 1)
457 cout <<
"starting tangent matrix computation" << endl;
458 pb.compute_tangent_matrix();
459 if (iter.get_noisy() > 1)
460 cout <<
"starting linear solver" << endl;
461 pb.linear_solve(dr, iter_linsolv);
462 if (!iter_linsolv.converged()) {
464 if (is_singular <= 4) {
465 if (iter.get_noisy())
466 cout <<
"Singular tangent matrix:"
467 " perturbation of the state vector." << endl;
469 pb.compute_residual();
471 if (iter.get_noisy())
472 cout <<
"Singular tangent matrix: perturbation failed, "
473 <<
"aborting." << endl;
477 else is_singular = 0;
479 if (iter.get_noisy() > 1) cout <<
"linear solver done" << endl;
482 pb.compute_residual();
484 R dec = R(1)/R(2), coeff2 = coeff * R(1.5);
486 while (dec > R(1E-5) && res >= res0 * coeff) {
487 gmm::add(gmm::scaled(dr, -dec), pb.state_vector());
488 pb.compute_residual();
490 if (res2 < res*R(0.95) || res2 >= res0 * coeff2) {
491 dec /= R(2); res = res2; coeff2 *= R(1.5);
493 gmm::add(gmm::scaled(dr, dec), pb.state_vector());
500 coeff = std::max(R(1.05), coeff*R(0.93));
501 bool near_end = (iter.get_iteration() > iter.get_maxiter()/2);
502 bool cut = (
alpha < R(1)) && near_end;
503 if ((res > minres && nit > 4) || cut) {
506 if ((stagn > 10 && alpha > alpha0 + R(5E-2)) || cut) {
508 if (near_end)
alpha = R(1);
510 pb.compute_residual();
513 stagn = 0; coeff = R(2);
516 if (res < minres || (alpha == R(1) && nit < 10)) stagn = 0;
518 if (nit < 5) minres = res0;
else minres = std::min(minres, res0);
520 if (iter.get_noisy())
521 cout <<
"step control [" << std::setw(8) << alpha0 <<
","
522 << std::setw(8) <<
alpha <<
"," << std::setw(10) << dec <<
"]";
526 crit = res / approx_eln;
529 if (iter.finished(crit)) {
530 if (iter.converged() && alpha < R(1)) {
532 alpha = std::min(R(1), alpha*R(3) - alpha0*R(2));
535 nit = 0; stagn = 0; coeff = R(2);
547 template <
typename PB>
550 typedef typename gmm::linalg_traits<typename PB::VECTOR>::value_type T;
551 typedef typename gmm::number_traits<T>::magnitude_type R;
553 iter_linsolv0.reduce_noisy();
554 iter_linsolv0.set_resmax(iter.get_resmax()/20.0);
555 iter_linsolv0.set_maxiter(10000);
557 pb.compute_residual();
558 R approx_eln = pb.approx_external_load_norm();
559 if (approx_eln == R(0)) approx_eln = R(1);
561 typename PB::VECTOR dr(gmm::vect_size(pb.residual()));
563 scalar_type crit = pb.residual_norm() / approx_eln;
564 while (!iter.finished(crit)) {
568 while (is_singular) {
571 if (iter.get_noisy() > 1)
572 cout <<
"starting computing tangent matrix" << endl;
573 pb.compute_tangent_matrix();
574 if (iter.get_noisy() > 1)
575 cout <<
"starting linear solver" << endl;
576 pb.linear_solve(dr, iter_linsolv);
577 if (!iter_linsolv.converged()) {
579 if (is_singular <= 4) {
580 if (iter.get_noisy())
581 cout <<
"Singular tangent matrix:"
582 " perturbation of the state vector." << endl;
584 pb.compute_residual();
586 if (iter.get_noisy())
587 cout <<
"Singular tangent matrix: perturbation failed, aborting."
592 else is_singular = 0;
595 if (iter.get_noisy() > 1) cout <<
"linear solver done" << endl;
596 R
alpha = pb.line_search(dr, iter);
598 if (iter.get_noisy()) cout <<
"alpha = " << std::setw(6) <<
alpha <<
" ";
600 crit = std::min(pb.residual_norm() / approx_eln,
601 gmm::vect_norm1(dr) / std::max(1E-25, pb.state_norm()));
611 typedef std::shared_ptr<abstract_linear_solver<model_real_sparse_matrix,
612 model_real_plain_vector> >
613 rmodel_plsolver_type;
614 typedef std::shared_ptr<abstract_linear_solver<model_complex_sparse_matrix,
615 model_complex_plain_vector> >
616 cmodel_plsolver_type;
618 template<
typename MATRIX,
typename VECTOR>
619 std::shared_ptr<abstract_linear_solver<MATRIX, VECTOR>>
620 default_linear_solver(
const model &md) {
622 #if GETFEM_PARA_LEVEL == 1 && GETFEM_PARA_SOLVER == MUMPS_PARA_SOLVER
623 if (md.is_symmetric())
624 return std::make_shared<linear_solver_mumps_sym<MATRIX, VECTOR>>();
626 return std::make_shared<linear_solver_mumps<MATRIX, VECTOR>>();
627 #elif GETFEM_PARA_LEVEL > 1 && GETFEM_PARA_SOLVER == MUMPS_PARA_SOLVER
628 if (md.is_symmetric())
629 return std::make_shared
630 <linear_solver_distributed_mumps_sym<MATRIX, VECTOR>>();
632 return std::make_shared
633 <linear_solver_distributed_mumps<MATRIX, VECTOR>>();
635 size_type ndof = md.nb_dof(), max3d = 15000, dim = md.leading_dimension();
636 # if defined(GMM_USES_MUMPS)
639 if ((ndof<300000 && dim<=2) || (ndof<max3d && dim<=3) || (ndof<1000)) {
640 # if defined(GMM_USES_MUMPS)
641 if (md.is_symmetric())
642 return std::make_shared<linear_solver_mumps_sym<MATRIX, VECTOR>>();
644 return std::make_shared<linear_solver_mumps<MATRIX, VECTOR>>();
645 # elif defined(GMM_USES_SUPERLU)
646 return std::make_shared<linear_solver_superlu<MATRIX, VECTOR>>();
649 "At least one direct solver (MUMPS or SuperLU) is required");
653 if (md.is_coercive())
654 return std::make_shared
655 <linear_solver_cg_preconditioned_ildlt<MATRIX, VECTOR>>();
658 return std::make_shared
659 <linear_solver_gmres_preconditioned_ilut<MATRIX,VECTOR>>();
661 return std::make_shared
662 <linear_solver_gmres_preconditioned_ilu<MATRIX,VECTOR>>();
666 return std::shared_ptr<abstract_linear_solver<MATRIX, VECTOR>>();
669 template <
typename MATRIX,
typename VECTOR>
670 std::shared_ptr<abstract_linear_solver<MATRIX, VECTOR>>
671 select_linear_solver(
const model &md,
const std::string &name) {
672 std::shared_ptr<abstract_linear_solver<MATRIX, VECTOR>> p;
673 if (bgeot::casecmp(name,
"superlu") == 0) {
674 #if defined(GMM_USES_SUPERLU)
675 return std::make_shared<linear_solver_superlu<MATRIX, VECTOR>>();
677 GMM_ASSERT1(
false,
"SuperLU is not interfaced");
680 else if (bgeot::casecmp(name,
"dense_lu") == 0)
681 return std::make_shared<linear_solver_dense_lu<MATRIX, VECTOR>>();
682 else if (bgeot::casecmp(name,
"mumps") == 0) {
683 #if defined(GMM_USES_MUMPS)
684 # if GETFEM_PARA_LEVEL <= 1
685 return std::make_shared<linear_solver_mumps<MATRIX, VECTOR>>();
687 return std::make_shared
688 <linear_solver_distributed_mumps<MATRIX, VECTOR>>();
691 GMM_ASSERT1(
false,
"Mumps is not interfaced");
694 else if (bgeot::casecmp(name,
"cg/ildlt") == 0)
695 return std::make_shared
696 <linear_solver_cg_preconditioned_ildlt<MATRIX, VECTOR>>();
697 else if (bgeot::casecmp(name,
"gmres/ilu") == 0)
698 return std::make_shared
699 <linear_solver_gmres_preconditioned_ilu<MATRIX, VECTOR>>();
700 else if (bgeot::casecmp(name,
"gmres/ilut") == 0)
701 return std::make_shared
702 <linear_solver_gmres_preconditioned_ilut<MATRIX, VECTOR>>();
703 else if (bgeot::casecmp(name,
"gmres/ilutp") == 0)
704 return std::make_shared
705 <linear_solver_gmres_preconditioned_ilutp<MATRIX, VECTOR>>();
706 else if (bgeot::casecmp(name,
"auto") == 0)
707 return default_linear_solver<MATRIX, VECTOR>(md);
709 GMM_ASSERT1(
false,
"Unknown linear solver");
710 return std::shared_ptr<abstract_linear_solver<MATRIX, VECTOR>>();
713 inline rmodel_plsolver_type rselect_linear_solver(
const model &md,
714 const std::string &name) {
715 return select_linear_solver<model_real_sparse_matrix,
716 model_real_plain_vector>(md, name);
719 inline cmodel_plsolver_type cselect_linear_solver(
const model &md,
720 const std::string &name) {
721 return select_linear_solver<model_complex_sparse_matrix,
722 model_complex_plain_vector>(md, name);
753 rmodel_plsolver_type lsolver,
754 abstract_newton_line_search &ls);
757 cmodel_plsolver_type lsolver,
758 abstract_newton_line_search &ls);
761 rmodel_plsolver_type lsolver);
764 cmodel_plsolver_type lsolver);
Incomplete Level 0 LDLT Preconditioner.
Incomplete LU without fill-in Preconditioner.
Incomplete LU with threshold and K fill-in Preconditioner.
ILUTP: Incomplete LU with threshold and K fill-in Preconditioner and column pivoting.
The Iteration object calculates whether the solution has reached the desired accuracy,...
Model representation in Getfem.
Interface with MUMPS (LU direct solver for sparse matrices).
void copy(const L1 &l1, L2 &l2)
*/
number_traits< typename linalg_traits< V1 >::value_type >::magnitude_type vect_dist1(const V1 &v1, const V2 &v2)
1-distance between two vectors
void clear(L &l)
clear (fill with zeros) a vector or matrix.
void add(const L1 &l1, L2 &l2)
*/
void lu_solve(const DenseMatrix &LU, const Pvector &pvector, VectorX &x, const VectorB &b)
LU Solve : Solve equation Ax=b, given an LU factored matrix.
Include standard gmm iterative solvers (cg, gmres, ...)
void gmres(const Mat &A, Vec &x, const VecB &b, const Precond &M, int restart, iteration &outer, Basis &KS)
Generalized Minimum Residual.
Interface with SuperLU (LU direct solver for sparse matrices).
size_t size_type
used as the common size type in the library
size_type alpha(short_type n, short_type d)
Return the value of which is the number of monomials of a polynomial of variables and degree .
GEneric Tool for Finite Element Methods.
void standard_solve(model &md, gmm::iteration &iter, rmodel_plsolver_type lsolver, abstract_newton_line_search &ls)
A default solver for the model brick system.