GetFEM  5.5
getfem_model_solvers.h
Go to the documentation of this file.
1 /* -*- c++ -*- (enables emacs c++ mode) */
2 /*===========================================================================
3 
4  Copyright (C) 2004-2026 Yves Renard
5 
6  This file is a part of GetFEM
7 
8  GetFEM is free software; you can redistribute it and/or modify it
9  under the terms of the GNU Lesser General Public License as published
10  by the Free Software Foundation; either version 3 of the License, or
11  (at your option) any later version along with the GCC Runtime Library
12  Exception either version 3.1 or (at your option) any later version.
13  This program is distributed in the hope that it will be useful, but
14  WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY
15  or FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
16  License and GCC Runtime Library Exception for more details.
17  You should have received a copy of the GNU Lesser General Public License
18  along with this program. If not, see https://www.gnu.org/licenses/.
19 
20  As a special exception, you may use this file as it is a part of a free
21  software library without restriction. Specifically, if other files
22  instantiate templates or use macros or inline functions from this file,
23  or you compile this file and link it with other files to produce an
24  executable, this file does not by itself cause the resulting executable
25  to be covered by the GNU Lesser General Public License. This exception
26  does not however invalidate any other reasons why the executable file
27  might be covered by the GNU Lesser General Public License.
28 
29 ===========================================================================*/
30 
31 /**
32  @file getfem_model_solvers.h
33  @author Yves Renard <Yves.Renard@insa-lyon.fr>
34  @date June 15, 2004.
35  @brief Standard solvers for model bricks
36  @see getfem_modeling.h
37 */
38 
39 #ifndef GETFEM_MODEL_SOLVERS_H__
40 #define GETFEM_MODEL_SOLVERS_H__
41 #include "getfem_models.h"
44 #include "gmm/gmm_iter.h"
45 #include "gmm/gmm_iter_solvers.h"
46 #include "gmm/gmm_dense_qr.h"
47 
48 //#include "gmm/gmm_inoutput.h"
49 
50 
51 namespace getfem {
52 
53  template <typename T> struct sort_abs_val_
54  { bool operator()(T x, T y) { return (gmm::abs(x) < gmm::abs(y)); } };
55 
56  template <typename MAT> void print_eigval(const MAT &M) {
57  // just to test a stiffness matrix if needing
58  typedef typename gmm::linalg_traits<MAT>::value_type T;
59  size_type n = gmm::mat_nrows(M);
60  gmm::dense_matrix<T> MM(n, n), Q(n, n);
61  std::vector<T> eigval(n);
62  gmm::copy(M, MM);
63  // gmm::symmetric_qr_algorithm(MM, eigval, Q);
64  gmm::implicit_qr_algorithm(MM, eigval, Q);
65  std::sort(eigval.begin(), eigval.end(), sort_abs_val_<T>());
66  cout << "eival = " << eigval << endl;
67 // cout << "vectp : " << gmm::mat_col(Q, n-1) << endl;
68 // cout << "vectp : " << gmm::mat_col(Q, n-2) << endl;
69 // double emax, emin;
70 // cout << "condition number" << condition_number(MM,emax,emin) << endl;
71 // cout << "emin = " << emin << endl;
72 // cout << "emax = " << emax << endl;
73  }
74 
75 
76  /* ***************************************************************** */
77  /* Linear solvers definition */
78  /* ***************************************************************** */
79 
80  template <typename MAT, typename VECT>
81  struct abstract_linear_solver {
82  typedef MAT MATRIX;
83  typedef VECT VECTOR;
84  virtual void operator ()(const MAT &, VECT &, const VECT &,
85  gmm::iteration &) const = 0;
86  virtual ~abstract_linear_solver() {}
87  };
88 
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,
93  gmm::iteration &iter) const {
95  gmm::cg(M, x, b, P, iter);
96  if (!iter.converged()) GMM_WARNING2("cg did not converge!");
97  }
98  };
99 
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,
104  gmm::iteration &iter) const {
106  gmm::gmres(M, x, b, P, 500, iter);
107  if (!iter.converged()) GMM_WARNING2("gmres did not converge!");
108  }
109  };
110 
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,
115  gmm::iteration &iter) const {
116  gmm::identity_matrix P;
117  gmm::gmres(M, x, b, P, 500, iter);
118  if (!iter.converged()) GMM_WARNING2("gmres did not converge!");
119  }
120  };
121 
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,
126  gmm::iteration &iter) const {
127  gmm::ilut_precond<MAT> P(M, 40, 1E-7);
128  gmm::gmres(M, x, b, P, 500, iter);
129  if (!iter.converged()) GMM_WARNING2("gmres did not converge!");
130  }
131  };
132 
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,
137  gmm::iteration &iter) const {
138  gmm::ilutp_precond<MAT> P(M, 20, 1E-7);
139  gmm::gmres(M, x, b, P, 500, iter);
140  if (!iter.converged()) GMM_WARNING2("gmres did not converge!");
141  }
142  };
143 
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,
149  gmm::iteration &iter) const {
150  double rcond;
151  /*gmm::HarwellBoeing_IO::write("test.hb", M);
152  std::fstream f("bbb", std::ios::out);
153  for (unsigned i=0; i < gmm::vect_size(b); ++i) f << b[i] << "\n";*/
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;
157  }
158  };
159 #endif
160 
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,
164  gmm::iteration &iter) const {
165  typedef typename gmm::linalg_traits<MAT>::value_type T;
166  gmm::dense_matrix<T> MM(gmm::mat_nrows(M),gmm::mat_ncols(M));
167  gmm::copy(M, MM);
168  gmm::lu_solve(MM, x, b);
169  iter.enforce_converged(true);
170  }
171  };
172 
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,
177  gmm::iteration &iter) const {
178  bool ok = gmm::MUMPS_solve(M, x, b, false);
179  iter.enforce_converged(ok);
180  }
181  };
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,
185  gmm::iteration &iter) const {
186  bool ok = gmm::MUMPS_solve(M, x, b, true);
187  iter.enforce_converged(ok);
188  }
189  };
190 #endif
191 
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,
197  gmm::iteration &iter) const {
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;
202  }
203  };
204 
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,
209  gmm::iteration &iter) const {
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;
214  }
215  };
216 #endif
217 
218 
219  /* ***************************************************************** */
220  /* Newton Line search definition */
221  /* ***************************************************************** */
222 
223  struct abstract_newton_line_search {
224  double conv_alpha, conv_r;
225  size_t it, itmax, glob_it;
226  // size_t tot_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) {
231  // tot_it += it; cout << "tot_it = " << tot_it << endl; it = 0;
232  return conv_alpha;
233  };
234  virtual double converged_residual(void) { return conv_r; };
235  virtual ~abstract_newton_line_search() { }
236  };
237 
238  // Dummy linesearch for the newton with step control
239  struct newton_search_with_step_control : public abstract_newton_line_search {
240 
241  virtual void init_search(double /*r*/, size_t /*git*/, double /*R0*/ = 0.0)
242  { GMM_ASSERT1(false, "Not to be used"); }
243 
244  virtual double next_try(void)
245  { GMM_ASSERT1(false, "Not to be used"); }
246 
247  virtual bool is_converged(double /*r*/, double /*R1*/ = 0.0)
248  { GMM_ASSERT1(false, "Not to be used"); }
249 
250  newton_search_with_step_control() {}
251  };
252 
253 
254  struct quadratic_newton_line_search : public abstract_newton_line_search {
255  double R0_, R1_;
256 
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");
260  glob_it = git;
261  conv_alpha = alpha = double(1); conv_r = first_res = r; it = 0;
262  R0_ = R0;
263  }
264  virtual double next_try(void) {
265  ++it;
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;
270  }
271  virtual bool is_converged(double r, double R1 = 0.0) {
272  conv_r = r;
273  R1_ = R1; return ((gmm::abs(R1_) < gmm::abs(R0_*0.5)) || it >= itmax);
274  }
275  quadratic_newton_line_search(size_t imax = size_t(-1)) { itmax = imax; }
276  };
277 
278 
279  struct simplest_newton_line_search : public abstract_newton_line_search {
280  double alpha, alpha_mult, first_res, alpha_max_ratio, alpha_min,
281  alpha_threshold_res;
282  virtual void init_search(double r, size_t git, double = 0.0) {
283  glob_it = git;
284  conv_alpha = alpha = double(1); conv_r = first_res = r; it = 0;
285  }
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) {
289  conv_r = r;
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)
293  || it >= itmax);
294  }
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)
301  { itmax = imax; }
302  };
303 
304  struct default_newton_line_search : public abstract_newton_line_search {
305  // This line search try to detect where is the minimum, dividing the step
306  // by a factor two each time.
307  // - it stops if the residual is less than the previous residual
308  // times alpha_min_ratio (= 0.9).
309  // - if the minimal step is reached with a residual greater than
310  // the previous residual times alpha_min_ratio then it decides
311  // between two options :
312  // - return with the step corresponding to the smallest residual
313  // - return with a greater residual corresponding to a residual
314  // less than the previous residual times alpha_max_ratio.
315  // the decision is taken regarding the previous iterations.
316  // - in order to shorten the line search, the process stops when
317  // the residual increases three times consecutively.
318  // possible improvment : detect the curvature at the origin
319  // (only one evaluation) and take it into account.
320  // Fitted to some experiments in contrib/tests_newton
321 
322  double alpha, alpha_old, alpha_mult, first_res, alpha_max_ratio;
323  double alpha_min_ratio, alpha_min;
324  size_type count, count_pat;
325  bool max_ratio_reached;
326  double alpha_max_ratio_reached, r_max_ratio_reached;
327  size_type it_max_ratio_reached;
328 
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; }
333  };
334 
335  /* the former default_newton_line_search */
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) {
340  glob_it = git;
341  conv_alpha = alpha = double(1);
342  prev_res = conv_r = first_res = r; it = 0;
343  }
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)
349  || it >= itmax)
350  { conv_r = r; return true; }
351  if (it > 1 && r > prev_res && prev_res < alpha_max_ratio * first_res)
352  return true;
353  prev_res = conv_r = r;
354  return false;
355  }
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; }
362  };
363 
364 
365  struct systematic_newton_line_search : public abstract_newton_line_search {
366  double alpha, alpha_mult, first_res;
367  double alpha_min, prev_res;
368  bool first;
369  virtual void init_search(double r, size_t git, double = 0.0) {
370  glob_it = git;
371  conv_alpha = alpha = double(1);
372  prev_res = conv_r = first_res = r; it = 0; first = true;
373  }
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) {
377  // cout << "a = " << alpha / alpha_mult << " r = " << r << endl;
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;
381  return false;
382  }
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; }
387  };
388 
389 
390  /* ***************************************************************** */
391  /* Newton(-Raphson) algorithm with step control. */
392  /* ***************************************************************** */
393  /* Still solves a problem F(X) = 0 starting at X0 but setting */
394  /* B0 = F(X0) and solving */
395  /* F(X) = (1-alpha)B0 (1) */
396  /* with alpha growing from 0 to 1. */
397  /* A very basic line search is applied. */
398  /* */
399  /* Step 0 : set alpha0 = 0, alpha = 1, X0 given and B0 = F(X0). */
400  /* Step 1 : Set Ri = F(Xi) - (1-alpha)B0 */
401  /* If ||Ri|| < rho, Xi+1 = Xi and go to step 2 */
402  /* Perform Newton step on problem (1) */
403  /* If the Newton do not converge (stagnation) */
404  /* alpha <- max(alpha0+1E-4, (alpha+alpha0)/2) */
405  /* Loop on step 1 with Xi unchanged */
406  /* Step 2 : if alpha=1 stop */
407  /* else alpha0 <- alpha, */
408  /* alpha <- min(1,alpha0+2*(alpha-alpha0)), */
409  /* Go to step 1 with Xi+1 */
410  /* being the result of Newton iterations of step1. */
411  /* */
412  /*********************************************************************/
413 
414  template <typename PB>
415  void Newton_with_step_control(PB &pb, gmm::iteration &iter)
416  {
417  typedef typename gmm::linalg_traits<typename PB::VECTOR>::value_type T;
418  typedef typename gmm::number_traits<T>::magnitude_type R;
419  gmm::iteration iter_linsolv0 = iter;
420  iter_linsolv0.reduce_noisy();
421  iter_linsolv0.set_resmax(iter.get_resmax()/20.0);
422  iter_linsolv0.set_maxiter(10000); // arbitrary
423 
424  pb.compute_residual();
425  R approx_eln = pb.approx_external_load_norm();
426  if (approx_eln == R(0)) approx_eln = R(1);
427 
428  typename PB::VECTOR b0(gmm::vect_size(pb.residual()));
429  gmm::copy(pb.residual(), b0);
430  typename PB::VECTOR Xi(gmm::vect_size(pb.residual()));
431  gmm::copy(pb.state_vector(), Xi);
432 
433  typename PB::VECTOR dr(gmm::vect_size(pb.residual()));
434 
435  R alpha0(0), alpha(1), res0(gmm::vect_norm1(b0)), minres(res0);
436  // const newton_search_with_step_control *ls
437  // = dynamic_cast<const newton_search_with_step_control *>(&(pb.ls));
438  // GMM_ASSERT1(ls, "Internal error");
439  size_type nit = 0, stagn = 0;
440  R coeff = R(2);
441 
442  scalar_type crit = pb.residual_norm() / approx_eln;
443  if (iter.finished(crit)) return;
444  for(;;) {
445 
446  crit = gmm::vect_dist1(pb.residual(), gmm::scaled(b0, R(1)-alpha))
447  / approx_eln;
448  if (!iter.converged(crit)) {
449  gmm::iteration iter_linsolv = iter_linsolv0;
450 
451  int is_singular = 1;
452  while (is_singular) { // Linear system solve
453  gmm::clear(dr);
454  pb.add_to_residual(b0, alpha-R(1)); // canceled at next compute_residual
455  iter_linsolv.init();
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()) {
463  is_singular++;
464  if (is_singular <= 4) {
465  if (iter.get_noisy())
466  cout << "Singular tangent matrix:"
467  " perturbation of the state vector." << endl;
468  pb.perturbation();
469  pb.compute_residual(); // cancels add_to_residual
470  } else {
471  if (iter.get_noisy())
472  cout << "Singular tangent matrix: perturbation failed, "
473  << "aborting." << endl;
474  return;
475  }
476  }
477  else is_singular = 0;
478  }
479  if (iter.get_noisy() > 1) cout << "linear solver done" << endl;
480 
481  gmm::add(dr, pb.state_vector());
482  pb.compute_residual(); // cancels add_to_residual
483  R res = gmm::vect_dist1(pb.residual(), gmm::scaled(b0, R(1)-alpha));
484  R dec = R(1)/R(2), coeff2 = coeff * R(1.5);
485 
486  while (dec > R(1E-5) && res >= res0 * coeff) {
487  gmm::add(gmm::scaled(dr, -dec), pb.state_vector());
488  pb.compute_residual();
489  R res2 = gmm::vect_dist1(pb.residual(), gmm::scaled(b0, R(1)-alpha));
490  if (res2 < res*R(0.95) || res2 >= res0 * coeff2) {
491  dec /= R(2); res = res2; coeff2 *= R(1.5);
492  } else {
493  gmm::add(gmm::scaled(dr, dec), pb.state_vector());
494  break;
495  }
496  }
497  dec *= R(2);
498 
499  nit++;
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) {
504  stagn++;
505 
506  if ((stagn > 10 && alpha > alpha0 + R(5E-2)) || cut) {
507  alpha = (alpha + alpha0) / R(2);
508  if (near_end) alpha = R(1);
509  gmm::copy(Xi, pb.state_vector());
510  pb.compute_residual();
511  res = gmm::vect_dist1(pb.residual(), gmm::scaled(b0, R(1)-alpha));
512  nit = 0;
513  stagn = 0; coeff = R(2);
514  }
515  }
516  if (res < minres || (alpha == R(1) && nit < 10)) stagn = 0;
517  res0 = res;
518  if (nit < 5) minres = res0; else minres = std::min(minres, res0);
519 
520  if (iter.get_noisy())
521  cout << "step control [" << std::setw(8) << alpha0 << ","
522  << std::setw(8) << alpha << "," << std::setw(10) << dec << "]";
523  ++iter;
524  // crit = std::min(res / approx_eln,
525  // gmm::vect_norm1(dr) / std::max(1E-25,pb.state_norm()));
526  crit = res / approx_eln;
527  }
528 
529  if (iter.finished(crit)) {
530  if (iter.converged() && alpha < R(1)) {
531  R a = alpha;
532  alpha = std::min(R(1), alpha*R(3) - alpha0*R(2));
533  alpha0 = a;
534  gmm::copy(pb.state_vector(), Xi);
535  nit = 0; stagn = 0; coeff = R(2);
536  } else return;
537  }
538  }
539  }
540 
541 
542 
543  /* ***************************************************************** */
544  /* Classical Newton(-Raphson) algorithm. */
545  /* ***************************************************************** */
546 
547  template <typename PB>
548  void classical_Newton(PB &pb, gmm::iteration &iter)
549  {
550  typedef typename gmm::linalg_traits<typename PB::VECTOR>::value_type T;
551  typedef typename gmm::number_traits<T>::magnitude_type R;
552  gmm::iteration iter_linsolv0 = iter;
553  iter_linsolv0.reduce_noisy();
554  iter_linsolv0.set_resmax(iter.get_resmax()/20.0);
555  iter_linsolv0.set_maxiter(10000); // arbitrary
556 
557  pb.compute_residual();
558  R approx_eln = pb.approx_external_load_norm();
559  if (approx_eln == R(0)) approx_eln = R(1);
560 
561  typename PB::VECTOR dr(gmm::vect_size(pb.residual()));
562 
563  scalar_type crit = pb.residual_norm() / approx_eln;
564  while (!iter.finished(crit)) {
565  gmm::iteration iter_linsolv = iter_linsolv0;
566 
567  int is_singular = 1;
568  while (is_singular) {
569  gmm::clear(dr);
570  iter_linsolv.init();
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()) {
578  is_singular++;
579  if (is_singular <= 4) {
580  if (iter.get_noisy())
581  cout << "Singular tangent matrix:"
582  " perturbation of the state vector." << endl;
583  pb.perturbation();
584  pb.compute_residual();
585  } else {
586  if (iter.get_noisy())
587  cout << "Singular tangent matrix: perturbation failed, aborting."
588  << endl;
589  return;
590  }
591  }
592  else is_singular = 0;
593  }
594 
595  if (iter.get_noisy() > 1) cout << "linear solver done" << endl;
596  R alpha = pb.line_search(dr, iter); //it is assumed that the linesearch
597  //executes a pb.compute_residual();
598  if (iter.get_noisy()) cout << "alpha = " << std::setw(6) << alpha << " ";
599  ++iter;
600  crit = std::min(pb.residual_norm() / approx_eln,
601  gmm::vect_norm1(dr) / std::max(1E-25, pb.state_norm()));
602  }
603  }
604 
605 
606 
607  //---------------------------------------------------------------------
608  // Default linear solver.
609  //---------------------------------------------------------------------
610 
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;
617 
618  template<typename MATRIX, typename VECTOR>
619  std::shared_ptr<abstract_linear_solver<MATRIX, VECTOR>>
620  default_linear_solver(const model &md) {
621 
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>>();
625  else
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>>();
631  else
632  return std::make_shared
633  <linear_solver_distributed_mumps<MATRIX, VECTOR>>();
634 #else
635  size_type ndof = md.nb_dof(), max3d = 15000, dim = md.leading_dimension();
636 # if defined(GMM_USES_MUMPS)
637  max3d = 250000;
638 # endif
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>>();
643  else
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>>();
647 # else
648  static_assert(false,
649  "At least one direct solver (MUMPS or SuperLU) is required");
650 # endif
651  }
652  else {
653  if (md.is_coercive())
654  return std::make_shared
655  <linear_solver_cg_preconditioned_ildlt<MATRIX, VECTOR>>();
656  else {
657  if (dim <= 2)
658  return std::make_shared
659  <linear_solver_gmres_preconditioned_ilut<MATRIX,VECTOR>>();
660  else
661  return std::make_shared
662  <linear_solver_gmres_preconditioned_ilu<MATRIX,VECTOR>>();
663  }
664  }
665 #endif
666  return std::shared_ptr<abstract_linear_solver<MATRIX, VECTOR>>();
667  }
668 
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>>();
676 #else
677  GMM_ASSERT1(false, "SuperLU is not interfaced");
678 #endif
679  }
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>>();
686 # else
687  return std::make_shared
688  <linear_solver_distributed_mumps<MATRIX, VECTOR>>();
689 # endif
690 #else
691  GMM_ASSERT1(false, "Mumps is not interfaced");
692 #endif
693  }
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);
708  else
709  GMM_ASSERT1(false, "Unknown linear solver");
710  return std::shared_ptr<abstract_linear_solver<MATRIX, VECTOR>>();
711  }
712 
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);
717  }
718 
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);
723  }
724 
725  //---------------------------------------------------------------------
726  // Standard solve.
727  //---------------------------------------------------------------------
728 
729 
730  /** A default solver for the model brick system.
731  Of course it could be not very well suited for a particular
732  problem, so it could be copied and adapted to change solvers,
733  add a special traitement on the problem, etc ... This is in
734  fact a model for your own solver.
735 
736  For small problems, a direct solver is used (gmm::SuperLU_solve),
737  for larger problems, a conjugate gradient gmm::cg (if the problem is
738  coercive) or a gmm::gmres is used (preconditioned with an incomplete
739  factorization).
740 
741  When MPI/METIS is enabled, a partition is done via METIS, and a parallel
742  solver can be used.
743 
744  Note that it is possible to disable some variables
745  (see the md.disable_variable(varname) method) in order to
746  solve the problem only with respect to a subset of variables (the
747  disabled variables are the considered as data) for instance to
748  replace the global Newton strategy with a fixed point one.
749 
750  @ingroup bricks
751  */
752  void standard_solve(model &md, gmm::iteration &iter,
753  rmodel_plsolver_type lsolver,
754  abstract_newton_line_search &ls);
755 
756  void standard_solve(model &md, gmm::iteration &iter,
757  cmodel_plsolver_type lsolver,
758  abstract_newton_line_search &ls);
759 
760  void standard_solve(model &md, gmm::iteration &iter,
761  rmodel_plsolver_type lsolver);
762 
763  void standard_solve(model &md, gmm::iteration &iter,
764  cmodel_plsolver_type lsolver);
765 
766  void standard_solve(model &md, gmm::iteration &iter);
767 
768 } /* end of namespace getfem. */
769 
770 
771 #endif /* GETFEM_MODEL_SOLVERS_H__ */
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,...
Definition: gmm_iter.h:52
Model representation in Getfem.
Interface with MUMPS (LU direct solver for sparse matrices).
void copy(const L1 &l1, L2 &l2)
*‍/
Definition: gmm_blas.h:976
number_traits< typename linalg_traits< V1 >::value_type >::magnitude_type vect_dist1(const V1 &v1, const V2 &v2)
1-distance between two vectors
Definition: gmm_blas.h:657
void clear(L &l)
clear (fill with zeros) a vector or matrix.
Definition: gmm_blas.h:58
void add(const L1 &l1, L2 &l2)
*‍/
Definition: gmm_blas.h:1275
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.
Definition: gmm_dense_lu.h:129
Dense QR factorization.
Iteration object.
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
Definition: bgeot_poly.h:48
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 .
Definition: bgeot_poly.cc:46
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.