38 #ifndef GMM_SOLVERS_SCHWARZ_ADDITIVE_H__ 
   39 #define GMM_SOLVERS_SCHWARZ_ADDITIVE_H__ 
   42 #if defined(GMM_USES_SUPERLU) 
   59   struct using_gmres {};
 
   60   struct using_bicgstab {};
 
   63   template <
typename P, 
typename local_solver, 
typename Matrix>
 
   64   struct actual_precond {
 
   66     static const APrecond &transform(
const P &PP) { 
return PP; }
 
   69   template <
typename Matrix1, 
typename Precond, 
typename Vector>
 
   70   void AS_local_solve(using_cg, 
const Matrix1 &A, Vector &x, 
const Vector &b,
 
   71                       const Precond &P, iteration &iter)
 
   72   { cg(A, x, b, P, iter); }
 
   74   template <
typename Matrix1, 
typename Precond, 
typename Vector>
 
   75   void AS_local_solve(using_gmres, 
const Matrix1 &A, Vector &x,
 
   76                       const Vector &b, 
const Precond &P, iteration &iter)
 
   77   { 
gmres(A, x, b, P, 100, iter); }
 
   79   template <
typename Matrix1, 
typename Precond, 
typename Vector>
 
   80   void AS_local_solve(using_bicgstab, 
const Matrix1 &A, Vector &x,
 
   81                       const Vector &b, 
const Precond &P, iteration &iter)
 
   82   { bicgstab(A, x, b, P, iter); }
 
   84   template <
typename Matrix1, 
typename Precond, 
typename Vector>
 
   85   void AS_local_solve(using_qmr, 
const Matrix1 &A, Vector &x,
 
   86                       const Vector &b, 
const Precond &P, iteration &iter)
 
   87   { 
qmr(A, x, b, P, iter); }
 
   89 #if defined(GMM_USES_SUPERLU) 
   90   struct using_superlu {};
 
   92   template <
typename P, 
typename Matrix>
 
   93   struct actual_precond<P, using_superlu, Matrix> {
 
   94     typedef typename linalg_traits<Matrix>::value_type value_type;
 
   95     typedef SuperLU_factor<value_type> APrecond;
 
   96     template <
typename PR>
 
   97     static APrecond transform(
const PR &) { 
return APrecond(); }
 
   98     static const APrecond &transform(
const APrecond &PP) { 
return PP; }
 
  101   template <
typename Matrix1, 
typename Precond, 
typename Vector>
 
  102   void AS_local_solve(using_superlu, 
const Matrix1 &, Vector &x,
 
  103                       const Vector &b, 
const Precond &P, iteration &iter)
 
  104   { P.solve(x, b); iter.set_iteration(1); }
 
  111   template <
typename Matrix1, 
typename Matrix2, 
typename Precond,
 
  112             typename local_solver>
 
  113   struct add_schwarz_mat{
 
  114     typedef typename linalg_traits<Matrix1>::value_type value_type;
 
  117     const std::vector<Matrix2> *vB;
 
  118     std::vector<Matrix2> vAloc;
 
  119     mutable iteration iter;
 
  122     mutable std::vector<std::vector<value_type> > gi, fi;
 
  123     std::vector<
typename actual_precond<Precond, local_solver,
 
  124                                         Matrix1>::APrecond> precond1;
 
  126     void init(
const Matrix1 &A_, 
const std::vector<Matrix2> &vB_,
 
  127               iteration iter_, 
const Precond &P, 
double residual_);
 
  130     add_schwarz_mat(
const Matrix1 &A_, 
const std::vector<Matrix2> &vB_,
 
  131                     iteration iter_, 
const Precond &P, 
double residual_)
 
  132     { init(A_, vB_, iter_, P, residual_); }
 
  135   template <
typename Matrix1, 
typename Matrix2, 
typename Precond,
 
  136             typename local_solver>
 
  137   void add_schwarz_mat<Matrix1, Matrix2, Precond, local_solver>::init(
 
  138        const Matrix1 &A_, 
const std::vector<Matrix2> &vB_,
 
  139        iteration iter_, 
const Precond &P, 
double residual_) {
 
  141     vB = &vB_; A = &A_; iter = iter_;
 
  142     residual = residual_;
 
  145     vAloc.resize(nb_sub);
 
  146     gi.resize(nb_sub); fi.resize(nb_sub);
 
  147     precond1.resize(nb_sub);
 
  148     std::fill(precond1.begin(), precond1.end(),
 
  149               actual_precond<Precond, local_solver, Matrix1>::transform(P));
 
  152     if (iter.get_noisy()) cout << 
"Init pour sub dom ";
 
  154     int size,tranche,borne_sup,borne_inf,rank,tag1=11,tag2=12,tag3=13,sizepr = 0;
 
  156     double t_ref,t_final;
 
  159     MPI_Comm_rank(MPI_COMM_WORLD, &rank);
 
  160     MPI_Comm_size(MPI_COMM_WORLD, &size);
 
  162     borne_inf=rank*tranche;
 
  163     borne_sup=(rank+1)*tranche;
 
  166     cout << 
"Nombre de sous domaines " << borne_sup - borne_inf << endl;
 
  168     int sizeA = mat_nrows(*A);
 
  169     gmm::csr_matrix<value_type> Acsr(sizeA, sizeA), Acsrtemp(sizeA, sizeA);
 
  171     int next = (rank + 1) % size;
 
  172     int previous = (rank + size - 1) % size;
 
  176     for (
int nproc = 0; nproc < size; ++nproc) {
 
  182         cout << 
"Sous domaines " << i << 
" : " << mat_ncols((*vB)[i]) << endl;
 
  186           if (iter.get_noisy())
 
  187             cout << i << 
" " << std::flush;
 
  188           Matrix2 Maux(mat_ncols((*vB)[i]), mat_nrows((*vB)[i]));
 
  191           Matrix2 Maux2(mat_ncols((*vB)[i]), mat_ncols((*vB)[i]));
 
  193             gmm::resize(vAloc[i], mat_ncols((*vB)[i]), mat_ncols((*vB)[i]));
 
  196           gmm::mult(gmm::transposed((*vB)[i]), Acsr, Maux);
 
  200           gmm::resize(vAloc[i], mat_ncols((*vB)[i]), mat_ncols((*vB)[i]));
 
  201           gmm::mult(gmm::transposed((*vB)[i]), *A, Maux);
 
  206           if (nproc == size - 1 ) {
 
  208             precond1[i].build_with(vAloc[i]);
 
  218       if (nproc != size - 1) {
 
  219         MPI_Sendrecv(&(Acsr.jc[0]), sizeA+1, MPI_INT, next, tag2,
 
  220                      &(Acsrtemp.jc[0]), sizeA+1, MPI_INT, previous, tag2,
 
  221                      MPI_COMM_WORLD, &status);
 
  222         if (Acsrtemp.jc[sizeA] > 
size_type(sizepr)) {
 
  223           sizepr = Acsrtemp.jc[sizeA];
 
  227         MPI_Sendrecv(&(Acsr.ir[0]), Acsr.jc[sizeA], MPI_INT, next, tag1,
 
  228                      &(Acsrtemp.ir[0]), Acsrtemp.jc[sizeA], MPI_INT, previous, tag1,
 
  229                      MPI_COMM_WORLD, &status);
 
  231         MPI_Sendrecv(&(Acsr.pr[0]), Acsr.jc[sizeA], mpi_type(value_type()), next, tag3,
 
  232                      &(Acsrtemp.pr[0]), Acsrtemp.jc[sizeA], mpi_type(value_type()), previous, tag3,
 
  233                      MPI_COMM_WORLD, &status);
 
  238     cout<<
"temps boucle precond "<< t_final-t_ref<<endl;
 
  240     if (iter.get_noisy()) cout << 
"\n";
 
  243   template <
typename Matrix1, 
typename Matrix2, 
typename Precond,
 
  244             typename Vector2, 
typename Vector3, 
typename local_solver>
 
  245   void mult(
const add_schwarz_mat<Matrix1, Matrix2, Precond, local_solver> &M,
 
  246             const Vector2 &p, Vector3 &q) {
 
  249     static double tmult_tot = 0.0;
 
  250     double t_ref = MPI_Wtime();
 
  255     tmult_tot += MPI_Wtime()-t_ref;
 
  256     cout << 
"tmult_tot = " << tmult_tot << endl;
 
  258     std::vector<double> qbis(gmm::vect_size(q));
 
  259     std::vector<double> qter(gmm::vect_size(q));
 
  264     int size,tranche,borne_sup,borne_inf,rank;
 
  266     MPI_Comm_rank(MPI_COMM_WORLD, &rank);
 
  267     MPI_Comm_size(MPI_COMM_WORLD, &size);
 
  269     borne_inf=rank*tranche;
 
  270     borne_sup=(rank+1)*tranche;
 
  279     for (
size_type i = 0; i < M.fi.size(); ++i)
 
  285         gmm::mult(gmm::transposed((*(M.vB))[i]), q, M.fi[i]);
 
  287        AS_local_solve(local_solver(), (M.vAloc)[i], (M.gi)[i],
 
  288                       (M.fi)[i],(M.precond1)[i],M.iter);
 
  289        itebilan = std::max(itebilan, M.iter.get_iteration());
 
  293     cout << 
"First  AS loop time " <<  MPI_Wtime() - t_ref << endl;
 
  303       for (
size_type i = 0; i < M.gi.size(); ++i)
 
  321     cout << 
"Second AS loop time " <<  MPI_Wtime() - t_ref << endl;
 
  327     static double t_tot = 0.0;
 
  346     MPI_Allreduce(&(qbis[0]), &(q[0]),gmm::vect_size(q), MPI_DOUBLE,
 
  347                   MPI_SUM,MPI_COMM_WORLD);
 
  349     t_tot += t_final-t_ref;
 
  350      cout<<
"["<< rank<<
"] temps reduce Resol "<< t_final-t_ref << 
" t_tot = " << t_tot << endl;
 
  353     if (M.iter.get_noisy() > 0) cout << 
"itebloc = " << itebilan << endl;
 
  354     M.itebilan += itebilan;
 
  355     M.iter.set_resmax((M.iter.get_resmax() + M.residual) * 0.5);
 
  358   template <
typename Matrix1, 
typename Matrix2, 
typename Precond,
 
  359             typename Vector2, 
typename Vector3, 
typename local_solver>
 
  360   void mult(
const add_schwarz_mat<Matrix1, Matrix2, Precond, local_solver> &M,
 
  361             const Vector2 &p, 
const Vector3 &q) {
 
  362     mult(M, p, 
const_cast<Vector3 &
>(q));
 
  365   template <
typename Matrix1, 
typename Matrix2, 
typename Precond,
 
  366             typename Vector2, 
typename Vector3, 
typename Vector4,
 
  367             typename local_solver>
 
  368   void mult(
const add_schwarz_mat<Matrix1, Matrix2, Precond, local_solver> &M,
 
  369             const Vector2 &p, 
const Vector3 &p2, Vector4 &q)
 
  372   template <
typename Matrix1, 
typename Matrix2, 
typename Precond,
 
  373             typename Vector2, 
typename Vector3, 
typename Vector4,
 
  374             typename local_solver>
 
  375   void mult(
const add_schwarz_mat<Matrix1, Matrix2, Precond, local_solver> &M,
 
  376             const Vector2 &p, 
const Vector3 &p2, 
const Vector4 &q)
 
  377   { 
mult(M, p, 
const_cast<Vector4 &
>(q)); 
add(p2, q); }
 
  383   template <
typename ASM_type, 
typename Vect>
 
  384   void AS_global_solve(using_cg, 
const ASM_type &ASM, Vect &x,
 
  385                        const Vect &b, iteration &iter)
 
  386   { cg(ASM, x, b, *(ASM.A), identity_matrix(), iter); }
 
  388   template <
typename ASM_type, 
typename Vect>
 
  389   void AS_global_solve(using_gmres, 
const ASM_type &ASM, Vect &x,
 
  390                        const Vect &b, iteration &iter)
 
  391   { 
gmres(ASM, x, b, identity_matrix(), 100, iter); }
 
  393   template <
typename ASM_type, 
typename Vect>
 
  394   void AS_global_solve(using_bicgstab, 
const ASM_type &ASM, Vect &x,
 
  395                        const Vect &b, iteration &iter)
 
  396   { bicgstab(ASM, x, b, identity_matrix(), iter); }
 
  398   template <
typename ASM_type, 
typename Vect>
 
  399   void AS_global_solve(using_qmr,
const ASM_type &ASM, Vect &x,
 
  400                        const Vect &b, iteration &iter)
 
  401   { 
qmr(ASM, x, b, identity_matrix(), iter); }
 
  403 #if defined(GMM_USES_SUPERLU) 
  404   template <
typename ASM_type, 
typename Vect>
 
  405   void AS_global_solve(using_superlu, 
const ASM_type &, Vect &,
 
  406                        const Vect &, iteration &) {
 
  407     GMM_ASSERT1(
false, 
"You cannot use SuperLU as " 
  408                 "global solver in additive Schwarz meethod");
 
  423   template <
typename Matrix1, 
typename Matrix2,
 
  424             typename Vector2, 
typename Vector3, 
typename Precond,
 
  425             typename local_solver, 
typename global_solver>
 
  427     add_schwarz_mat<Matrix1, Matrix2, Precond, local_solver> &ASM, Vector3 &u,
 
  428     const Vector2 &f, 
iteration &iter, 
const global_solver&) {
 
  430     typedef typename linalg_traits<Matrix1>::value_type value_type;
 
  432     size_type nb_sub = ASM.vB->size(), nb_dof = gmm::vect_size(f);
 
  434     std::vector<value_type> g(nb_dof);
 
  435     std::vector<value_type> gbis(nb_dof);
 
  437     double t_init=MPI_Wtime();
 
  438     int size,tranche,borne_sup,borne_inf,rank;
 
  439     MPI_Comm_rank(MPI_COMM_WORLD, &rank);
 
  440     MPI_Comm_size(MPI_COMM_WORLD, &size);
 
  442     borne_inf=rank*tranche;
 
  443     borne_sup=(rank+1)*tranche;
 
  450     for (size_type i = 0; i < nb_sub; ++i)
 
  457       gmm::mult(gmm::transposed((*(ASM.vB))[i]), f, ASM.fi[i]);
 
  459       AS_local_solve(local_solver(), ASM.vAloc[i], ASM.gi[i], ASM.fi[i],
 
  460                      ASM.precond1[i], ASM.iter);
 
  461       ASM.itebilan = std::max(ASM.itebilan, ASM.iter.get_iteration());
 
  463     gmm::mult((*(ASM.vB))[i], ASM.gi[i], gbis,gbis);
 
  465     gmm::mult((*(ASM.vB))[i], ASM.gi[i], g, g);
 
  469     cout<<
"temps boucle init "<< MPI_Wtime()-t_init<<endl;
 
  470     double t_ref,t_final;
 
  472     MPI_Allreduce(&(gbis[0]), &(g[0]),gmm::vect_size(g), MPI_DOUBLE,
 
  473                   MPI_SUM,MPI_COMM_WORLD);
 
  475     cout<<
"temps reduce init "<< t_final-t_ref<<endl;
 
  479     cout<<
"begin global AS"<<endl;
 
  481     AS_global_solve(global_solver(), ASM, u, g, iter);
 
  484     cout<<
"temps AS Global Solve "<< t_final-t_ref<<endl;
 
  486     if (iter.get_noisy())
 
  487       cout << 
"Total number of internal iterations : " << ASM.itebilan << endl;
 
  493   template <
typename Matrix1, 
typename Matrix2,
 
  494             typename Vector2, 
typename Vector3, 
typename Precond,
 
  495             typename local_solver, 
typename global_solver>
 
  497                         const Vector2 &f, 
const Precond &P,
 
  498                         const std::vector<Matrix2> &vB,
 
  500                         local_solver, global_solver) {
 
  502     if (iter.get_rhsnorm() == 0.0) { 
gmm::clear(u); 
return; }
 
  503     iteration iter2 = iter; iter2.reduce_noisy();
 
  505     add_schwarz_mat<Matrix1, Matrix2, Precond, local_solver>
 
  506       ASM(A, vB, iter2, P, iter.get_resmax());
 
  518   template <
typename Matrixt, 
typename MatrixBi>
 
  519   class NewtonAS_struct {
 
  522     typedef Matrixt tangent_matrix_type;
 
  523     typedef MatrixBi B_matrix_type;
 
  524     typedef typename linalg_traits<Matrixt>::value_type value_type;
 
  525     typedef std::vector<value_type> Vector;
 
  528     virtual const std::vector<MatrixBi> &get_vB() = 0;
 
  530     virtual void compute_F(Vector &f, Vector &x) = 0;
 
  531     virtual void compute_tangent_matrix(Matrixt &M, Vector &x) = 0;
 
  533     virtual void compute_sub_tangent_matrix(Matrixt &Mloc, Vector &x,
 
  536     virtual void compute_sub_F(Vector &fi, Vector &x, 
size_type i) = 0;
 
  538     virtual ~NewtonAS_struct() {}
 
  541   template <
typename Matrixt, 
typename MatrixBi>
 
  542   struct AS_exact_gradient {
 
  543     const std::vector<MatrixBi> &vB;
 
  544     std::vector<Matrixt> vM;
 
  545     std::vector<Matrixt> vMloc;
 
  548       for (
size_type i = 0; i < vB.size(); ++i) {
 
  549         Matrixt aux(gmm::mat_ncols(vB[i]), gmm::mat_ncols(vM[i]));
 
  550         gmm::resize(vMloc[i], gmm::mat_ncols(vB[i]), gmm::mat_ncols(vB[i]));
 
  551         gmm::mult(gmm::transposed(vB[i]), vM[i], aux);
 
  555     AS_exact_gradient(
const std::vector<MatrixBi> &vB_) : vB(vB_) {
 
  556       vM.resize(vB.size()); vMloc.resize(vB.size());
 
  557       for (
size_type i = 0; i < vB.size(); ++i) {
 
  558         gmm::resize(vM[i], gmm::mat_nrows(vB[i]), gmm::mat_nrows(vB[i]));
 
  563   template <
typename Matrixt, 
typename MatrixBi,
 
  564             typename Vector2, 
typename Vector3>
 
  565   void mult(
const AS_exact_gradient<Matrixt, MatrixBi> &M,
 
  566             const Vector2 &p, Vector3 &q) {
 
  568     typedef typename gmm::linalg_traits<Vector3>::value_type T;
 
  569     std::vector<T> v(gmm::vect_size(p)), w, x;
 
  570     for (
size_type i = 0; i < M.vB.size(); ++i) {
 
  571       w.resize(gmm::mat_ncols(M.vB[i]));
 
  572       x.resize(gmm::mat_ncols(M.vB[i]));
 
  574       gmm::mult(gmm::transposed(M.vB[i]), v, w);
 
  575 #if defined(GMM_USES_SUPERLU) 
  577       gmm::SuperLU_solve(M.vMloc[i], x, w, rcond);
 
  579       gmm::MUMPS_solve(M.vMloc[i], x, w);
 
  587   template <
typename Matrixt, 
typename MatrixBi,
 
  588             typename Vector2, 
typename Vector3>
 
  589   void mult(
const AS_exact_gradient<Matrixt, MatrixBi> &M,
 
  590             const Vector2 &p, 
const Vector3 &q) {
 
  591     mult(M, p, 
const_cast<Vector3 &
>(q));
 
  594   template <
typename Matrixt, 
typename MatrixBi,
 
  595             typename Vector2, 
typename Vector3, 
typename Vector4>
 
  596   void mult(
const AS_exact_gradient<Matrixt, MatrixBi> &M,
 
  597             const Vector2 &p, 
const Vector3 &p2, Vector4 &q)
 
  600   template <
typename Matrixt, 
typename MatrixBi,
 
  601             typename Vector2, 
typename Vector3, 
typename Vector4>
 
  602   void mult(
const AS_exact_gradient<Matrixt, MatrixBi> &M,
 
  603             const Vector2 &p, 
const Vector3 &p2, 
const Vector4 &q)
 
  604   { 
mult(M, p, 
const_cast<Vector4 &
>(q)); 
add(p2, q); }
 
  606   struct S_default_newton_line_search {
 
  608     double conv_alpha, conv_r;
 
  609     size_t it, itmax, glob_it;
 
  611     double alpha, alpha_old, alpha_mult, first_res, alpha_max_ratio;
 
  612     double alpha_min_ratio, alpha_min;
 
  614     bool max_ratio_reached;
 
  615     double alpha_max_ratio_reached, r_max_ratio_reached;
 
  619     double converged_value() { 
return conv_alpha; };
 
  620     double converged_residual() { 
return conv_r; };
 
  622     virtual void init_search(
double r, 
size_t git, 
double = 0.0) {
 
  623       alpha_min_ratio = 0.9;
 
  625       alpha_max_ratio = 10.0;
 
  628       glob_it = git; 
if (git <= 1) count_pat = 0;
 
  629       conv_alpha = 
alpha = alpha_old = 1.;
 
  630       conv_r = first_res = r; it = 0;
 
  632       max_ratio_reached = 
false;
 
  634     virtual double next_try() {
 
  636       if (alpha >= 0.4) 
alpha *= 0.5; 
else alpha *= alpha_mult; ++it;
 
  639     virtual bool is_converged(
double r, 
double = 0.0) {
 
  641       if (!max_ratio_reached && r < first_res * alpha_max_ratio) {
 
  642         alpha_max_ratio_reached = alpha_old; r_max_ratio_reached = r;
 
  643         it_max_ratio_reached = it; max_ratio_reached = 
true;
 
  645       if (max_ratio_reached && r < r_max_ratio_reached * 0.5
 
  646           && r > first_res * 1.1 && it <= it_max_ratio_reached+1) {
 
  647         alpha_max_ratio_reached = alpha_old; r_max_ratio_reached = r;
 
  648         it_max_ratio_reached = it;
 
  650       if (count == 0 || r < conv_r)
 
  651         { conv_r = r; conv_alpha = alpha_old; count = 1; }
 
  652       if (conv_r < first_res) ++count;
 
  654       if (r < first_res *  alpha_min_ratio)
 
  655         { count_pat = 0; 
return true; }
 
  656       if (count >= 5 || (alpha < alpha_min && max_ratio_reached)) {
 
  657         if (conv_r < first_res * 0.99) count_pat = 0;
 
  659           { conv_r=r_max_ratio_reached; conv_alpha=alpha_max_ratio_reached; }
 
  660         if (conv_r >= first_res * 0.9999) count_pat++;
 
  665     S_default_newton_line_search() { count_pat = 0; }
 
  670   template <
typename Matrixt, 
typename MatrixBi, 
typename Vector,
 
  671             typename Precond, 
typename local_solver, 
typename global_solver>
 
  672   void Newton_additive_Schwarz(NewtonAS_struct<Matrixt, MatrixBi> &NS,
 
  674                                iteration &iter, 
const Precond &P,
 
  675                                local_solver, global_solver) {
 
  676     Vector &u = 
const_cast<Vector &
>(u_);
 
  677     typedef typename linalg_traits<Vector>::value_type value_type;
 
  678     typedef typename number_traits<value_type>::magnitude_type mtype;
 
  679     typedef actual_precond<Precond, local_solver, Matrixt> chgt_precond;
 
  681     double residual = iter.get_resmax();
 
  683     S_default_newton_line_search internal_ls;
 
  684     S_default_newton_line_search external_ls;
 
  686     typename chgt_precond::APrecond PP = chgt_precond::transform(P);
 
  687     iter.set_rhsnorm(mtype(1));
 
  688     iteration iternc(iter);
 
  689     iternc.reduce_noisy(); iternc.set_maxiter(
size_type(-1));
 
  690     iteration iter2(iternc);
 
  691     iteration iter3(iter2); iter3.reduce_noisy();
 
  692     iteration iter4(iter3);
 
  693     iternc.set_name(
"Local Newton");
 
  694     iter2.set_name(
"Linear System for Global Newton");
 
  695     iternc.set_resmax(residual/100.0);
 
  696     iter3.set_resmax(residual/10000.0);
 
  697     iter2.set_resmax(residual/1000.0);
 
  698     iter4.set_resmax(residual/1000.0);
 
  699     std::vector<value_type> rhs(NS.size()), x(NS.size()), d(NS.size());
 
  700     std::vector<value_type> xi, xii, fi, di;
 
  702     std::vector< std::vector<value_type> > vx(NS.get_vB().size());
 
  703     for (
size_type i = 0; i < NS.get_vB().size(); ++i) 
 
  706     Matrixt Mloc, M(NS.size(), NS.size());
 
  707     NS.compute_F(rhs, u);
 
  708     mtype act_res=
gmm::vect_norm2(rhs), act_res_new(0), precond_res = act_res;
 
  711     while(!iter.finished(std::min(act_res, precond_res))) {
 
  712       for (
int SOR_step = 0;  SOR_step >= 0; --SOR_step) {
 
  714         for (
size_type isd = 0; isd < NS.get_vB().size(); ++isd) {
 
  715           const MatrixBi &Bi = (NS.get_vB())[isd];
 
  718           xi.resize(si); xii.resize(si); fi.resize(si); di.resize(si);
 
  721           iternc.set_maxiter(30); 
 
  722           if (iternc.get_noisy())
 
  723             cout << 
"Non-linear local problem " << isd << endl;
 
  726           NS.compute_sub_F(fi, x, isd); gmm::scale(fi, value_type(-1));
 
  728           if (r > value_type(0)) {
 
  729             iternc.set_rhsnorm(std::max(r, mtype(1)));
 
  730             while(!iternc.finished(r)) {
 
  731               NS.compute_sub_tangent_matrix(Mloc, x, isd);
 
  735               AS_local_solve(local_solver(), Mloc, di, fi, PP, iter3);
 
  737               internal_ls.init_search(r, iternc.get_iteration());
 
  739                 alpha = internal_ls.next_try();
 
  740                 gmm::add(xi, gmm::scaled(di, -alpha), xii);
 
  741                 gmm::mult(Bi, gmm::scaled(xii, -1.0), u, x);
 
  742                 NS.compute_sub_F(fi, x, isd); gmm::scale(fi, value_type(-1));
 
  744               } 
while (!internal_ls.is_converged(r_t));
 
  746               if (alpha != internal_ls.converged_value()) {
 
  747                 alpha = internal_ls.converged_value();
 
  748                 gmm::add(xi, gmm::scaled(di, -alpha), xii);
 
  749                 gmm::mult(Bi, gmm::scaled(xii, -1.0), u, x);
 
  750                 NS.compute_sub_F(fi, x, isd); gmm::scale(fi, value_type(-1));
 
  755               if (iternc.get_noisy()) cout << 
"(step=" << 
alpha << 
")\t";
 
  758             if (SOR_step) 
gmm::mult(Bi, gmm::scaled(xii, -1.0), u, u);
 
  759             gmm::mult(Bi, gmm::scaled(xii, -1.0), rhs, rhs);
 
  763         if (SOR_step) cout << 
"SOR step residual = " << precond_res << endl;
 
  764         if (precond_res < residual) 
break;
 
  765         cout << 
"Precond residual = " << precond_res << endl;
 
  771         NS.compute_tangent_matrix(M, u);
 
  772         add_schwarz_mat<Matrixt, MatrixBi, Precond, local_solver>
 
  773           ASM(M, NS.get_vB(), iter4, P, iter.get_resmax());
 
  774         AS_global_solve(global_solver(), ASM, d, rhs, iter2);
 
  777         AS_exact_gradient<Matrixt, MatrixBi> eg(NS.get_vB());
 
  778         for (
size_type i = 0; i < NS.get_vB().size(); ++i) {
 
  779           NS.compute_tangent_matrix(eg.vM[i], vx[i]);
 
  782         gmres(eg, d, rhs, gmm::identity_matrix(), 50, iter2);
 
  786       external_ls.init_search(act_res, iter.get_iteration());
 
  788         alpha = external_ls.next_try();
 
  789         gmm::add(gmm::scaled(d, alpha), u, x);
 
  790         NS.compute_F(rhs, x);
 
  792       } 
while (!external_ls.is_converged(act_res_new));
 
  794       if (alpha != external_ls.converged_value()) {
 
  795         alpha = external_ls.converged_value();
 
  796         gmm::add(gmm::scaled(d, alpha), u, x);
 
  797         NS.compute_F(rhs, x);
 
  801       if (iter.get_noisy() > 1) cout << endl;
 
  802       act_res = act_res_new;
 
  803       if (iter.get_noisy())
 
  804         cout << 
"(step=" << 
alpha 
  805              << 
")\t unprecond res = " << act_res << 
" ";
 
The Iteration object calculates whether the solution has reached the desired accuracy,...
 
Interface with MUMPS (LU direct solver for sparse matrices).
 
void mult_add(const L1 &l1, const L2 &l2, L3 &l3)
*/
 
void copy(const L1 &l1, L2 &l2)
*/
 
number_traits< typename linalg_traits< V >::value_type >::magnitude_type vect_norm2(const V &v)
Euclidean norm of a vector.
 
void clear(L &l)
clear (fill with zeros) a vector or matrix.
 
void resize(V &v, size_type n)
*/
 
void mult(const L1 &l1, const L2 &l2, L3 &l3)
*/
 
void add(const L1 &l1, L2 &l2)
*/
 
Include the base gmm files.
 
void additive_schwarz(add_schwarz_mat< Matrix1, Matrix2, Precond, local_solver > &ASM, Vector3 &u, const Vector2 &f, iteration &iter, const global_solver &)
Function to call if the ASM matrix is precomputed for successive solve with the same system.
 
BiCGStab iterative solver.
 
Conjugate gradient iterative solver.
 
GMRES (Generalized Minimum Residual) iterative solver.
 
void gmres(const Mat &A, Vec &x, const VecB &b, const Precond &M, int restart, iteration &outer, Basis &KS)
Generalized Minimum Residual.
 
Quasi-Minimal Residual iterative solver.
 
void qmr(const Matrix &A, Vector &x, const VectorB &b, const Precond1 &M1, iteration &iter)
Quasi-Minimal 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 .