40 #ifndef GETFEM_MODEL_SOLVERS_H__ 
   41 #define GETFEM_MODEL_SOLVERS_H__ 
   54   template <
typename T> 
struct sort_abs_val_
 
   55   { 
bool operator()(T x, T y) { 
return (gmm::abs(x) < gmm::abs(y)); } };
 
   57   template <
typename MAT> 
void print_eigval(
const MAT &M) {
 
   59     typedef typename gmm::linalg_traits<MAT>::value_type T;
 
   61     gmm::dense_matrix<T> MM(n, n), Q(n, n);
 
   62     std::vector<T> eigval(n);
 
   65     gmm::implicit_qr_algorithm(MM, eigval, Q);
 
   66     std::sort(eigval.begin(), eigval.end(), sort_abs_val_<T>());
 
   67     cout << 
"eival = " << eigval << endl;
 
   81   template <
typename MAT, 
typename VECT>
 
   82   struct abstract_linear_solver {
 
   85     virtual void operator ()(
const MAT &, VECT &, 
const VECT &,
 
   87     virtual ~abstract_linear_solver() {}
 
   90   template <
typename MAT, 
typename VECT>
 
   91   struct linear_solver_cg_preconditioned_ildlt
 
   92     : 
public abstract_linear_solver<MAT, VECT> {
 
   93     void operator ()(
const MAT &M, VECT &x, 
const VECT &b,
 
   96       gmm::cg(M, x, b, P, iter);
 
   97       if (!iter.converged()) GMM_WARNING2(
"cg did not converge!");
 
  101   template <
typename MAT, 
typename VECT>
 
  102   struct linear_solver_gmres_preconditioned_ilu
 
  103     : 
public abstract_linear_solver<MAT, VECT> {
 
  104     void operator ()(
const MAT &M, VECT &x, 
const VECT &b,
 
  108       if (!iter.converged()) GMM_WARNING2(
"gmres did not converge!");
 
  112   template <
typename MAT, 
typename VECT>
 
  113   struct linear_solver_gmres_unpreconditioned
 
  114     : 
public abstract_linear_solver<MAT, VECT> {
 
  115     void operator ()(
const MAT &M, VECT &x, 
const VECT &b,
 
  117       gmm::identity_matrix P;
 
  119       if (!iter.converged()) GMM_WARNING2(
"gmres did not converge!");
 
  123   template <
typename MAT, 
typename VECT>
 
  124   struct linear_solver_gmres_preconditioned_ilut
 
  125     : 
public abstract_linear_solver<MAT, VECT> {
 
  126     void operator ()(
const MAT &M, VECT &x, 
const VECT &b,
 
  130       if (!iter.converged()) GMM_WARNING2(
"gmres did not converge!");
 
  134   template <
typename MAT, 
typename VECT>
 
  135   struct linear_solver_gmres_preconditioned_ilutp
 
  136     : 
public abstract_linear_solver<MAT, VECT> {
 
  137     void operator ()(
const MAT &M, VECT &x, 
const VECT &b,
 
  141       if (!iter.converged()) GMM_WARNING2(
"gmres did not converge!");
 
  145 #if defined(GMM_USES_SUPERLU) 
  146   template <
typename MAT, 
typename VECT>
 
  147   struct linear_solver_superlu
 
  148     : 
public abstract_linear_solver<MAT, VECT> {
 
  149     void operator ()(
const MAT &M, VECT &x, 
const VECT &b,
 
  155       int info = gmm::SuperLU_solve(M, x, b, rcond);
 
  156       iter.enforce_converged(info == 0);
 
  157       if (iter.get_noisy()) cout << 
"condition number: " << 1.0/rcond<< endl;
 
  162   template <
typename MAT, 
typename VECT>
 
  163   struct linear_solver_dense_lu : 
public abstract_linear_solver<MAT, VECT> {
 
  164     void operator ()(
const MAT &M, VECT &x, 
const VECT &b,
 
  166       typedef typename gmm::linalg_traits<MAT>::value_type T;
 
  167       gmm::dense_matrix<T> MM(gmm::mat_nrows(M),gmm::mat_ncols(M));
 
  170       iter.enforce_converged(
true);
 
  174 #if defined(GMM_USES_MUMPS) 
  175   template <
typename MAT, 
typename VECT>
 
  176   struct linear_solver_mumps : 
public abstract_linear_solver<MAT, VECT> {
 
  177     void operator ()(
const MAT &M, VECT &x, 
const VECT &b,
 
  179       bool ok = gmm::MUMPS_solve(M, x, b, 
false);
 
  180       iter.enforce_converged(ok);
 
  183   template <
typename MAT, 
typename VECT>
 
  184   struct linear_solver_mumps_sym : 
public abstract_linear_solver<MAT, VECT> {
 
  185     void operator ()(
const MAT &M, VECT &x, 
const VECT &b,
 
  187       bool ok = gmm::MUMPS_solve(M, x, b, 
true);
 
  188       iter.enforce_converged(ok);
 
  193 #if GETFEM_PARA_LEVEL > 1 && GETFEM_PARA_SOLVER == MUMPS_PARA_SOLVER 
  194   template <
typename MAT, 
typename VECT>
 
  195   struct linear_solver_distributed_mumps
 
  196     : 
public abstract_linear_solver<MAT, VECT> {
 
  197     void operator ()(
const MAT &M, VECT &x, 
const VECT &b,
 
  199       double tt_ref=MPI_Wtime();
 
  200       bool ok = MUMPS_distributed_matrix_solve(M, x, b, 
false);
 
  201       iter.enforce_converged(ok);
 
  202       if (MPI_IS_MASTER()) cout<<
"UNSYMMETRIC MUMPS time "<< MPI_Wtime() - tt_ref<<endl;
 
  206   template <
typename MAT, 
typename VECT>
 
  207   struct linear_solver_distributed_mumps_sym
 
  208     : 
public abstract_linear_solver<MAT, VECT> {
 
  209     void operator ()(
const MAT &M, VECT &x, 
const VECT &b,
 
  211       double tt_ref=MPI_Wtime();
 
  212       bool ok = MUMPS_distributed_matrix_solve(M, x, b, 
true);
 
  213       iter.enforce_converged(ok);
 
  214       if (MPI_IS_MASTER()) cout<<
"SYMMETRIC MUMPS time "<< MPI_Wtime() - tt_ref<<endl;
 
  224   struct abstract_newton_line_search {
 
  225     double conv_alpha, conv_r;
 
  226     size_t it, itmax, glob_it;
 
  228     virtual void init_search(
double r, 
size_t git, 
double R0 = 0.0) = 0;
 
  229     virtual double next_try(
void) = 0;
 
  230     virtual bool is_converged(
double, 
double R1 = 0.0) = 0;
 
  231     virtual double converged_value(
void) {
 
  235     virtual double converged_residual(
void) { 
return conv_r; };
 
  236     virtual ~abstract_newton_line_search() { }
 
  240   struct newton_search_with_step_control : 
public abstract_newton_line_search {
 
  242     virtual void init_search(
double , 
size_t , 
double  = 0.0)
 
  243     { GMM_ASSERT1(
false, 
"Not to be used"); }
 
  245     virtual double next_try(
void)
 
  246     { GMM_ASSERT1(
false, 
"Not to be used"); }
 
  248     virtual bool is_converged(
double , 
double  = 0.0)
 
  249     { GMM_ASSERT1(
false, 
"Not to be used"); }
 
  251     newton_search_with_step_control() {}
 
  255   struct quadratic_newton_line_search : 
public abstract_newton_line_search {
 
  258     double alpha, first_res;
 
  259     virtual void init_search(
double r, 
size_t git, 
double R0 = 0.0) {
 
  260       GMM_ASSERT1(R0 != 0.0, 
"You have to specify R0");
 
  262       conv_alpha = 
alpha = double(1); conv_r = first_res = r; it = 0;
 
  265     virtual double next_try(
void) {
 
  267       if (it == 1) 
return double(1);
 
  268       GMM_ASSERT1(R1_ != 0.0, 
"You have to specify R1");
 
  269       double a = R0_ / R1_;
 
  270       return (a < 0) ? (a*0.5 + sqrt(a*a*0.25-a)) : a*0.5;
 
  272     virtual bool is_converged(
double r, 
double R1 = 0.0) {
 
  274       R1_ = R1; 
return ((gmm::abs(R1_) < gmm::abs(R0_*0.5)) || it >= itmax);
 
  276     quadratic_newton_line_search(
size_t imax = 
size_t(-1))  { itmax = imax; }
 
  280   struct simplest_newton_line_search : 
public abstract_newton_line_search {
 
  281     double alpha, alpha_mult, first_res, alpha_max_ratio, alpha_min,
 
  283     virtual void init_search(
double r, 
size_t git, 
double = 0.0) {
 
  285       conv_alpha = 
alpha = double(1); conv_r = first_res = r; it = 0;
 
  287     virtual double next_try(
void)
 
  288     { conv_alpha = 
alpha; 
alpha *= alpha_mult; ++it; 
return conv_alpha; }
 
  289     virtual bool is_converged(
double r, 
double = 0.0) {
 
  291       return ((it <= 1 && r < first_res)
 
  292               || (r <= first_res * alpha_max_ratio && r <= alpha_threshold_res)
 
  293               || (conv_alpha <= alpha_min && r < first_res * 1e5)
 
  296     simplest_newton_line_search
 
  297     (
size_t imax = 
size_t(-1), 
double a_max_ratio = 6.0/5.0,
 
  298      double a_min = 1.0/1000.0, 
double a_mult = 3.0/5.0,
 
  299      double a_threshold_res = 1e50)
 
  300       : alpha_mult(a_mult), alpha_max_ratio(a_max_ratio), alpha_min(a_min),
 
  301         alpha_threshold_res(a_threshold_res)
 
  305   struct default_newton_line_search : 
public abstract_newton_line_search {
 
  323     double alpha, alpha_old, alpha_mult, first_res, alpha_max_ratio;
 
  324     double alpha_min_ratio, alpha_min;
 
  326     bool max_ratio_reached;
 
  327     double alpha_max_ratio_reached, r_max_ratio_reached;
 
  330     virtual void init_search(
double r, 
size_t git, 
double = 0.0);
 
  331     virtual double next_try(
void);
 
  332     virtual bool is_converged(
double r, 
double = 0.0);
 
  333     default_newton_line_search(
void) { count_pat = 0; }
 
  337   struct basic_newton_line_search : 
public abstract_newton_line_search {
 
  338     double alpha, alpha_mult, first_res, alpha_max_ratio;
 
  339     double alpha_min, prev_res, alpha_max_augment;
 
  340     virtual void init_search(
double r, 
size_t git, 
double = 0.0) {
 
  342       conv_alpha = 
alpha = double(1);
 
  343       prev_res = conv_r = first_res = r; it = 0;
 
  345     virtual double next_try(
void)
 
  346     { conv_alpha = 
alpha; 
alpha *= alpha_mult; ++it; 
return conv_alpha; }
 
  347     virtual bool is_converged(
double r, 
double = 0.0) {
 
  348       if (glob_it == 0 || (r < first_res / 
double(2))
 
  349           || (conv_alpha <= alpha_min && r < first_res * alpha_max_augment)
 
  351         { conv_r = r; 
return true; }
 
  352       if (it > 1 && r > prev_res && prev_res < alpha_max_ratio * first_res)
 
  354       prev_res = conv_r = r;
 
  357     basic_newton_line_search
 
  358     (
size_t imax = 
size_t(-1),
 
  359      double a_max_ratio = 5.0/3.0,
 
  360      double a_min = 1.0/1000.0, 
double a_mult = 3.0/5.0, 
double a_augm = 2.0)
 
  361       : alpha_mult(a_mult), alpha_max_ratio(a_max_ratio),
 
  362         alpha_min(a_min), alpha_max_augment(a_augm) { itmax = imax; }
 
  366   struct systematic_newton_line_search : 
public abstract_newton_line_search {
 
  367     double alpha, alpha_mult, first_res;
 
  368     double alpha_min, prev_res;
 
  370     virtual void init_search(
double r, 
size_t git, 
double = 0.0) {
 
  372       conv_alpha = 
alpha = double(1);
 
  373       prev_res = conv_r = first_res = r; it = 0; first = 
true;
 
  375     virtual double next_try(
void)
 
  376     { 
double a = 
alpha; 
alpha *= alpha_mult; ++it; 
return a; }
 
  377     virtual bool is_converged(
double r, 
double = 0.0) {
 
  379       if (r < conv_r || first)
 
  380         { conv_r = r; conv_alpha = 
alpha / alpha_mult; first = 
false; }
 
  381       if ((alpha <= alpha_min*alpha_mult) || it >= itmax) 
return true;
 
  384     systematic_newton_line_search
 
  385     (
size_t imax = 
size_t(-1),
 
  386      double a_min = 1.0/10000.0, 
double a_mult = 3.0/5.0)
 
  387       : alpha_mult(a_mult), alpha_min(a_min)  { itmax = imax; }
 
  415   template <
typename PB>
 
  418     typedef typename gmm::linalg_traits<typename PB::VECTOR>::value_type T;
 
  419     typedef typename gmm::number_traits<T>::magnitude_type R;
 
  421     iter_linsolv0.reduce_noisy();
 
  422     iter_linsolv0.set_resmax(iter.get_resmax()/20.0);
 
  423     iter_linsolv0.set_maxiter(10000); 
 
  425     pb.compute_residual();
 
  426     R approx_eln = pb.approx_external_load_norm();
 
  427     if (approx_eln == R(0)) approx_eln = R(1);
 
  429     typename PB::VECTOR b0(gmm::vect_size(pb.residual()));
 
  431     typename PB::VECTOR Xi(gmm::vect_size(pb.residual()));
 
  434     typename PB::VECTOR dr(gmm::vect_size(pb.residual()));
 
  436     R alpha0(0), 
alpha(1), res0(gmm::vect_norm1(b0)), minres(res0);
 
  443     scalar_type crit = pb.residual_norm() / approx_eln;
 
  444     if (iter.finished(crit)) 
return;
 
  449       if (!iter.converged(crit)) {
 
  453         while (is_singular) { 
 
  455           pb.add_to_residual(b0, alpha-R(1)); 
 
  457           if (iter.get_noisy() > 1)
 
  458             cout << 
"starting tangent matrix computation" << endl;
 
  459           pb.compute_tangent_matrix();
 
  460           if (iter.get_noisy() > 1)
 
  461             cout << 
"starting linear solver" << endl;
 
  462           pb.linear_solve(dr, iter_linsolv);
 
  463           if (!iter_linsolv.converged()) {
 
  465             if (is_singular <= 4) {
 
  466               if (iter.get_noisy())
 
  467                 cout << 
"Singular tangent matrix:" 
  468                   " perturbation of the state vector." << endl;
 
  470               pb.compute_residual(); 
 
  472               if (iter.get_noisy())
 
  473                 cout << 
"Singular tangent matrix: perturbation failed, " 
  474                      << 
"aborting." << endl;
 
  478           else is_singular = 0;
 
  480         if (iter.get_noisy() > 1) cout << 
"linear solver done" << endl;
 
  483         pb.compute_residual(); 
 
  485         R dec = R(1)/R(2), coeff2 = coeff * R(1.5);
 
  487         while (dec > R(1E-5) && res >= res0 * coeff) {
 
  488           gmm::add(gmm::scaled(dr, -dec), pb.state_vector());
 
  489           pb.compute_residual();
 
  491           if (res2 < res*R(0.95) || res2 >= res0 * coeff2) {
 
  492             dec /= R(2); res = res2; coeff2 *= R(1.5);
 
  494             gmm::add(gmm::scaled(dr, dec), pb.state_vector());
 
  501         coeff = std::max(R(1.05), coeff*R(0.93));
 
  502         bool near_end = (iter.get_iteration() > iter.get_maxiter()/2);
 
  503         bool cut = (
alpha < R(1)) && near_end;
 
  504         if ((res > minres && nit > 4) || cut) {
 
  507           if ((stagn > 10 && alpha > alpha0 + R(5E-2)) || cut) {
 
  509             if (near_end) 
alpha = R(1);
 
  511             pb.compute_residual();
 
  514             stagn = 0; coeff = R(2);
 
  517         if (res < minres || (alpha == R(1) &&  nit < 10)) stagn = 0;
 
  519         if (nit < 5) minres = res0; 
else minres = std::min(minres, res0);
 
  521         if (iter.get_noisy())
 
  522           cout << 
"step control [" << std::setw(8) << alpha0 << 
"," 
  523                << std::setw(8) << 
alpha << 
"," << std::setw(10) << dec <<  
"]";
 
  527         crit = res / approx_eln;
 
  530       if (iter.finished(crit)) {
 
  531         if (iter.converged() && alpha < R(1)) {
 
  533           alpha = std::min(R(1), alpha*R(3) - alpha0*R(2));
 
  536           nit = 0; stagn = 0; coeff = R(2);
 
  548   template <
typename PB>
 
  551     typedef typename gmm::linalg_traits<typename PB::VECTOR>::value_type T;
 
  552     typedef typename gmm::number_traits<T>::magnitude_type R;
 
  554     iter_linsolv0.reduce_noisy();
 
  555     iter_linsolv0.set_resmax(iter.get_resmax()/20.0);
 
  556     iter_linsolv0.set_maxiter(10000); 
 
  558     pb.compute_residual();
 
  559     R approx_eln = pb.approx_external_load_norm();
 
  560     if (approx_eln == R(0)) approx_eln = R(1);
 
  562     typename PB::VECTOR dr(gmm::vect_size(pb.residual()));
 
  564     scalar_type crit = pb.residual_norm() / approx_eln;
 
  565     while (!iter.finished(crit)) {
 
  569       while (is_singular) {
 
  572         if (iter.get_noisy() > 1)
 
  573           cout << 
"starting computing tangent matrix" << endl;
 
  574         pb.compute_tangent_matrix();
 
  575         if (iter.get_noisy() > 1)
 
  576           cout << 
"starting linear solver" << endl;
 
  577         pb.linear_solve(dr, iter_linsolv);
 
  578         if (!iter_linsolv.converged()) {
 
  580           if (is_singular <= 4) {
 
  581             if (iter.get_noisy())
 
  582               cout << 
"Singular tangent matrix:" 
  583                 " perturbation of the state vector." << endl;
 
  585             pb.compute_residual();
 
  587             if (iter.get_noisy())
 
  588               cout << 
"Singular tangent matrix: perturbation failed, aborting." 
  593         else is_singular = 0;
 
  596       if (iter.get_noisy() > 1) cout << 
"linear solver done" << endl;
 
  597       R 
alpha = pb.line_search(dr, iter); 
 
  599       if (iter.get_noisy()) cout << 
"alpha = " << std::setw(6) << 
alpha << 
" ";
 
  601       crit = std::min(pb.residual_norm() / approx_eln,
 
  602                       gmm::vect_norm1(dr) / std::max(1E-25, pb.state_norm()));
 
  612   typedef std::shared_ptr<abstract_linear_solver<model_real_sparse_matrix,
 
  613                                                  model_real_plain_vector> >
 
  614     rmodel_plsolver_type;
 
  615   typedef std::shared_ptr<abstract_linear_solver<model_complex_sparse_matrix,
 
  616                                                  model_complex_plain_vector> >
 
  617     cmodel_plsolver_type;
 
  619   template<
typename MATRIX, 
typename VECTOR>
 
  620   std::shared_ptr<abstract_linear_solver<MATRIX, VECTOR>>
 
  621   default_linear_solver(
const model &md) {
 
  623 #if GETFEM_PARA_LEVEL == 1 && GETFEM_PARA_SOLVER == MUMPS_PARA_SOLVER 
  624     if (md.is_symmetric())
 
  625       return std::make_shared<linear_solver_mumps_sym<MATRIX, VECTOR>>();
 
  627       return std::make_shared<linear_solver_mumps<MATRIX, VECTOR>>();
 
  628 #elif GETFEM_PARA_LEVEL > 1 && GETFEM_PARA_SOLVER == MUMPS_PARA_SOLVER 
  629     if (md.is_symmetric())
 
  630       return std::make_shared
 
  631         <linear_solver_distributed_mumps_sym<MATRIX, VECTOR>>();
 
  633         return std::make_shared
 
  634           <linear_solver_distributed_mumps<MATRIX, VECTOR>>();
 
  636     size_type ndof = md.nb_dof(), max3d = 15000, dim = md.leading_dimension();
 
  637 # if defined(GMM_USES_MUMPS) 
  640     if ((ndof<300000 && dim<=2) || (ndof<max3d && dim<=3) || (ndof<1000)) {
 
  641 # if defined(GMM_USES_MUMPS) 
  642       if (md.is_symmetric())
 
  643         return std::make_shared<linear_solver_mumps_sym<MATRIX, VECTOR>>();
 
  645         return std::make_shared<linear_solver_mumps<MATRIX, VECTOR>>();
 
  646 # elif defined(GMM_USES_SUPERLU) 
  647       return std::make_shared<linear_solver_superlu<MATRIX, VECTOR>>();
 
  650                     "At least one direct solver (MUMPS or SuperLU) is required");
 
  654       if (md.is_coercive())
 
  655         return std::make_shared
 
  656           <linear_solver_cg_preconditioned_ildlt<MATRIX, VECTOR>>();
 
  659           return std::make_shared
 
  660             <linear_solver_gmres_preconditioned_ilut<MATRIX,VECTOR>>();
 
  662             return std::make_shared
 
  663               <linear_solver_gmres_preconditioned_ilu<MATRIX,VECTOR>>();
 
  667     return std::shared_ptr<abstract_linear_solver<MATRIX, VECTOR>>();
 
  670   template <
typename MATRIX, 
typename VECTOR>
 
  671   std::shared_ptr<abstract_linear_solver<MATRIX, VECTOR>>
 
  672   select_linear_solver(
const model &md, 
const std::string &name) {
 
  673     std::shared_ptr<abstract_linear_solver<MATRIX, VECTOR>> p;
 
  674     if (bgeot::casecmp(name, 
"superlu") == 0) {
 
  675 #if defined(GMM_USES_SUPERLU) 
  676       return std::make_shared<linear_solver_superlu<MATRIX, VECTOR>>();
 
  678       GMM_ASSERT1(
false, 
"SuperLU is not interfaced");
 
  681     else if (bgeot::casecmp(name, 
"dense_lu") == 0)
 
  682       return std::make_shared<linear_solver_dense_lu<MATRIX, VECTOR>>();
 
  683     else if (bgeot::casecmp(name, 
"mumps") == 0) {
 
  684 #if defined(GMM_USES_MUMPS) 
  685 # if GETFEM_PARA_LEVEL <= 1 
  686       return std::make_shared<linear_solver_mumps<MATRIX, VECTOR>>();
 
  688       return std::make_shared
 
  689         <linear_solver_distributed_mumps<MATRIX, VECTOR>>();
 
  692       GMM_ASSERT1(
false, 
"Mumps is not interfaced");
 
  695     else if (bgeot::casecmp(name, 
"cg/ildlt") == 0)
 
  696       return std::make_shared
 
  697         <linear_solver_cg_preconditioned_ildlt<MATRIX, VECTOR>>();
 
  698     else if (bgeot::casecmp(name, 
"gmres/ilu") == 0)
 
  699       return std::make_shared
 
  700         <linear_solver_gmres_preconditioned_ilu<MATRIX, VECTOR>>();
 
  701     else if (bgeot::casecmp(name, 
"gmres/ilut") == 0)
 
  702       return std::make_shared
 
  703         <linear_solver_gmres_preconditioned_ilut<MATRIX, VECTOR>>();
 
  704     else if (bgeot::casecmp(name, 
"gmres/ilutp") == 0)
 
  705       return std::make_shared
 
  706         <linear_solver_gmres_preconditioned_ilutp<MATRIX, VECTOR>>();
 
  707     else if (bgeot::casecmp(name, 
"auto") == 0)
 
  708       return default_linear_solver<MATRIX, VECTOR>(md);
 
  710       GMM_ASSERT1(
false, 
"Unknown linear solver");
 
  711     return std::shared_ptr<abstract_linear_solver<MATRIX, VECTOR>>();
 
  714   inline rmodel_plsolver_type rselect_linear_solver(
const model &md,
 
  715                                              const std::string &name) {
 
  716     return select_linear_solver<model_real_sparse_matrix,
 
  717                                 model_real_plain_vector>(md, name);
 
  720   inline cmodel_plsolver_type cselect_linear_solver(
const model &md,
 
  721                                              const std::string &name) {
 
  722     return select_linear_solver<model_complex_sparse_matrix,
 
  723                                 model_complex_plain_vector>(md, name);
 
  754                       rmodel_plsolver_type lsolver,
 
  755                       abstract_newton_line_search &ls);
 
  758                       cmodel_plsolver_type lsolver,
 
  759                       abstract_newton_line_search &ls);
 
  762                       rmodel_plsolver_type lsolver);
 
  765                       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.