Bridge++  Version 1.5.4
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
solver_BiCGStab_DS_L_Cmplx.cpp
Go to the documentation of this file.
1 
15 
16 #ifdef USE_FACTORY_AUTOREGISTER
17 namespace {
18  bool init = Solver_BiCGStab_DS_L_Cmplx::register_factory();
19 }
20 #endif
21 
22 const std::string Solver_BiCGStab_DS_L_Cmplx::class_name = "Solver_BiCGStab_DS_L_Cmplx";
23 
24 //====================================================================
26 {
27  const string str_vlevel = params.get_string("verbose_level");
28 
29  m_vl = vout.set_verbose_level(str_vlevel);
30 
31  //- fetch and check input parameters
32  int Niter, Nrestart;
33  double Stop_cond;
34  bool use_init_guess;
35  int N_L;
36  double Omega_tolerance;
37  double Tol_L;
38 
39  int err = 0;
40  err += params.fetch_int("maximum_number_of_iteration", Niter);
41  err += params.fetch_int("maximum_number_of_restart", Nrestart);
42  err += params.fetch_double("convergence_criterion_squared", Stop_cond);
43  err += params.fetch_bool("use_initial_guess", use_init_guess);
44  err += params.fetch_double("Omega_tolerance", Omega_tolerance);
45  err += params.fetch_int("number_of_orthonormal_vectors", N_L);
46  err += params.fetch_double("tolerance_for_DynamicSelection_of_L", Tol_L);
47 
48  if (err) {
49  vout.crucial(m_vl, "Error at %s: input parameter not found.\n", class_name.c_str());
50  exit(EXIT_FAILURE);
51  }
52 
53 
54  set_parameters(Niter, Nrestart, Stop_cond, use_init_guess);
55  set_parameters_BiCGStab_series(Omega_tolerance);
56  set_parameters_DS_L(N_L, Tol_L);
57 }
58 
59 
60 //====================================================================
61 void Solver_BiCGStab_DS_L_Cmplx::set_parameters(const int Niter, const int Nrestart, const double Stop_cond)
62 {
64 
65  //- print input parameters
66  vout.general(m_vl, "%s: input parameters\n", class_name.c_str());
67  vout.general(m_vl, " Niter = %d\n", Niter);
68  vout.general(m_vl, " Nrestart = %d\n", Nrestart);
69  vout.general(m_vl, " Stop_cond = %8.2e\n", Stop_cond);
70 
71  //- range check
72  int err = 0;
73  err += ParameterCheck::non_negative(Niter);
74  err += ParameterCheck::non_negative(Nrestart);
75  err += ParameterCheck::square_non_zero(Stop_cond);
76 
77  if (err) {
78  vout.crucial(m_vl, "Error at %s: parameter range check failed.\n", class_name.c_str());
79  exit(EXIT_FAILURE);
80  }
81 
82  //- store values
83  m_Niter = Niter;
84  m_Nrestart = Nrestart;
85  m_Stop_cond = Stop_cond;
86 }
87 
88 
89 //====================================================================
90 void Solver_BiCGStab_DS_L_Cmplx::set_parameters(const int Niter, const int Nrestart, const double Stop_cond, const bool use_init_guess)
91 {
93 
94  //- print input parameters
95  vout.general(m_vl, "%s: input parameters\n", class_name.c_str());
96  vout.general(m_vl, " Niter = %d\n", Niter);
97  vout.general(m_vl, " Nrestart = %d\n", Nrestart);
98  vout.general(m_vl, " Stop_cond = %8.2e\n", Stop_cond);
99  vout.general(m_vl, " use_init_guess = %s\n", use_init_guess ? "true" : "false");
100 
101  //- range check
102  int err = 0;
103  err += ParameterCheck::non_negative(Niter);
104  err += ParameterCheck::non_negative(Nrestart);
105  err += ParameterCheck::square_non_zero(Stop_cond);
106 
107  if (err) {
108  vout.crucial(m_vl, "Error at %s: parameter range check failed.\n", class_name.c_str());
109  exit(EXIT_FAILURE);
110  }
111 
112  //- store values
113  m_Niter = Niter;
114  m_Nrestart = Nrestart;
115  m_Stop_cond = Stop_cond;
116  m_use_init_guess = use_init_guess;
117 }
118 
119 
120 //====================================================================
122 {
124 
125  //- print input parameters
126  vout.general(m_vl, " Omega_tolerance = %8.2e\n", Omega_tolerance);
127 
128  //- range check
129  // NB. Omega_tolerance == 0.0 is allowed.
130 
131  //- store values
132  m_Omega_tolerance = Omega_tolerance;
133 }
134 
135 
136 //====================================================================
137 void Solver_BiCGStab_DS_L_Cmplx::set_parameters_DS_L(const int N_L, const double Tol_L)
138 {
139  //- print input parameters
140  vout.general(m_vl, " N_L = %d\n", N_L);
141  vout.general(m_vl, " Tol_L = %16.8e\n", Tol_L);
142 
143  //- range check
144  int err = 0;
145  err += ParameterCheck::non_negative(N_L);
146  err += ParameterCheck::square_non_zero(Tol_L);
147 
148  if (err) {
149  vout.crucial(m_vl, "Error at %s: parameter range check failed.\n", class_name.c_str());
150  exit(EXIT_FAILURE);
151  }
152 
153  //- store values
154  m_N_L = N_L;
155  m_Tol_L = Tol_L;
156 }
157 
158 
159 //====================================================================
161  int& Nconv, double& diff)
162 {
163  const double bnorm2 = b.norm2();
164  const int bsize = b.size();
165 
166  vout.paranoiac(m_vl, "%s: starts\n", class_name.c_str());
167  vout.paranoiac(m_vl, " norm of b = %16.8e\n", bnorm2);
168  vout.paranoiac(m_vl, " size of b = %d\n", bsize);
169 
170  bool is_converged = false;
171  int Nconv2 = 0;
172  double diff2 = 1.0; // superficial initialization
173  double rr;
174 
175  int Nconv_unit = 1;
176  // if (m_fopr->get_mode() == "DdagD" || m_fopr->get_mode() == "DDdag") {
177  // Nconv_unit = 2;
178  // }
179 
180  reset_field(b);
181 
182  if (m_use_init_guess) {
183  copy(m_s, xq); // s = xq;
184  } else {
185  copy(m_s, b); // s = b;
186  }
187  solve_init(b, rr);
188  Nconv2 += Nconv_unit;
189 
190  vout.detailed(m_vl, " iter: %8d %22.15e\n", Nconv2, rr / bnorm2);
191 
192 
193  for (int i_restart = 0; i_restart < m_Nrestart; i_restart++) {
194  for (int iter = 0; iter < m_Niter; iter++) {
195  if (rr / bnorm2 < m_Stop_cond) break;
196 
197  solve_step(rr);
198  Nconv2 += 2 * Nconv_unit * m_N_L_prev;
199 
200  vout.paranoiac(m_vl, " iter,N_L: %8d %8d\n", Nconv2, m_N_L_prev);
201  vout.detailed(m_vl, " iter: %8d %22.15e\n", Nconv2, rr / bnorm2);
202  }
203 
204  //- calculate true residual
205  m_fopr->mult(m_s, m_x); // s = m_fopr->mult(x);
206  axpy(m_s, -1.0, b); // s -= b;
207  diff2 = m_s.norm2();
208 
209  if (diff2 / bnorm2 < m_Stop_cond) {
210  vout.detailed(m_vl, "%s: converged.\n", class_name.c_str());
211  vout.detailed(m_vl, " iter(final): %8d %22.15e\n", Nconv2, diff2 / bnorm2);
212 
213  is_converged = true;
214 
215  m_Nrestart_count = i_restart;
216  m_Nconv_count = Nconv2;
217 
218  break;
219  } else {
220  //- restart with new approximate solution
221  copy(m_s, m_x); // s = x;
222  solve_init(b, rr);
223 
224  vout.detailed(m_vl, "%s: restarted.\n", class_name.c_str());
225  }
226  }
227 
228 
229  if (!is_converged) {
230  vout.crucial(m_vl, "Error at %s: not converged.\n", class_name.c_str());
231  vout.crucial(m_vl, " iter(final): %8d %22.15e\n", Nconv2, diff2 / bnorm2);
232  exit(EXIT_FAILURE);
233  }
234 
235 
236  copy(xq, m_x); // xq = x;
237 
238 #pragma omp barrier
239 #pragma omp master
240  {
241  diff = sqrt(diff2 / bnorm2);
242  Nconv = Nconv2;
243  }
244 #pragma omp barrier
245 }
246 
247 
248 //====================================================================
250 {
251 #pragma omp barrier
252 #pragma omp master
253  {
254  const int Nin = b.nin();
255  const int Nvol = b.nvol();
256  const int Nex = b.nex();
257 
258  if ((m_s.nin() != Nin) || (m_s.nvol() != Nvol) || (m_s.nex() != Nex)) {
259  m_s.reset(Nin, Nvol, Nex);
260  m_x.reset(Nin, Nvol, Nex);
261  m_r_init.reset(Nin, Nvol, Nex);
262  m_v.reset(Nin, Nvol, Nex);
263  }
264 
265  m_u.resize(m_N_L + 1);
266  m_r.resize(m_N_L + 1);
267 
268  for (int i = 0; i < m_N_L + 1; ++i) {
269  m_u[i].reset(Nin, Nvol, Nex);
270  m_r[i].reset(Nin, Nvol, Nex);
271  }
272  }
273 #pragma omp barrier
274 
275  vout.paranoiac(m_vl, "%s: field size reset.\n", class_name.c_str());
276 }
277 
278 
279 //====================================================================
281 {
282  copy(m_x, m_s); // x = s;
283 
284  for (int i = 0; i < m_N_L + 1; ++i) {
285  m_r[i].set(0.0); // r[i] = 0.0;
286  m_u[i].set(0.0); // u[i] = 0.0;
287  }
288 
289  // r[0] = b - A x_0
290  m_fopr->mult(m_v, m_s); // v = m_fopr->mult(s);
291  copy(m_r[0], b); // r[0] = b;
292  axpy(m_r[0], -1.0, m_v); // r[0] -= v;
293 
294  copy(m_r_init, m_r[0]); // r_init = r[0];
295  rr = m_r[0].norm2(); // rr = r[0] * r[0];
296 
297 #pragma omp barrier
298 #pragma omp master
299  {
300  m_rho_prev = cmplx(-1.0, 0.0);
301 
302  // NB. m_alpha_prev = 0.0 \neq 1.0
303  m_alpha_prev = cmplx(0.0, 0.0);
304 
305  m_N_L_prev = m_N_L;
306  }
307 #pragma omp barrier
308 }
309 
310 
311 //====================================================================
313 {
314  dcomplex rho_prev2 = m_rho_prev;
315  dcomplex alpha_prev2 = m_alpha_prev;
316 
317  int N_L_tmp = 0; // superficial initialization
318  dcomplex c_Rayleigh_prev = cmplx(0.0, 0.0);
319 
320  bool is_converged_L = false;
321 
322 
323  for (int j = 0; j < m_N_L_prev; ++j) {
324  if (!is_converged_L) {
325  dcomplex rho = dotc(m_r[j], m_r_init); // rho = r[j] * r_init;
326  rho = conj(rho);
327 
328  dcomplex beta = alpha_prev2 * (rho / rho_prev2);
329 
330  rho_prev2 = rho;
331 
332  for (int i = 0; i < j + 1; ++i) {
333  aypx(-beta, m_u[i], m_r[i]); // u[i] = - beta * u[i] + r[i];
334  }
335 
336  m_fopr->mult(m_u[j + 1], m_u[j]); // u[j+1] = m_fopr->mult(u[j]);
337 
338  dcomplex gamma = dotc(m_u[j + 1], m_r_init);
339  alpha_prev2 = rho_prev2 / conj(gamma);
340 
341  for (int i = 0; i < j + 1; ++i) {
342  axpy(m_r[i], -alpha_prev2, m_u[i + 1]); // r[i] -= alpha_prev2 * u[i+1];
343  }
344 
345  m_fopr->mult(m_r[j + 1], m_r[j]); // r[j+1] = m_fopr->mult(r[j]);
346 
347  axpy(m_x, alpha_prev2, m_u[0]); // x += alpha_prev2 * u[0];
348 
349 
350  //- calculate the Rayleigh quotient for N_L_tmp
351 
352  // dcomplex c_Rayleigh = (r[j] * r[j+1]) / (r[j] * r[j]);
353  double r_tmp = m_r[j].norm2(); // r_tmp = r[j] * r[j];
354  dcomplex c_Rayleigh = dotc(m_r[j], m_r[j + 1]) / r_tmp; // c_Rayleigh = r[j] * r[j+1] / r_tmp;
355 
356  dcomplex c_E = (c_Rayleigh - c_Rayleigh_prev) / c_Rayleigh;
357 
358  c_Rayleigh_prev = c_Rayleigh;
359 
360  N_L_tmp = j + 1;
361 
362  // vout.paranoiac(m_vl, "N_L_tmp,abs(c_E),m_Tol_L = %d %f %f\n",N_L_tmp,abs(c_E),m_Tol_L);
363 
364  if (abs(c_E) < m_Tol_L) {
365  // vout.paranoiac(m_vl, "N_L_tmp = %d\n",N_L_tmp);
366 
367  // break;
368  is_converged_L = true;
369  }
370  }
371  }
372 
373 
374  std::vector<double> sigma(m_N_L + 1);
375  std::vector<dcomplex> gamma_prime(m_N_L + 1);
376 
377  // NB. tau(m_N_L,m_N_L+1), not (m_N_L+1,m_N_L+1)
378  std::vector<dcomplex> tau(m_N_L * (m_N_L + 1));
379 
380  const double sigma_0 = m_r[0].norm2();
381 
382  for (int j = 1; j < N_L_tmp + 1; ++j) {
383  for (int i = 1; i < j; ++i) {
384  int ij = index_ij(i, j);
385 
386  dcomplex r_ji = dotc(m_r[j], m_r[i]);
387  tau[ij] = conj(r_ji) / sigma[i]; // tau[ij] = (r[j] * r[i]) / sigma[i];
388  axpy(m_r[j], -tau[ij], m_r[i]); // r[j] -= tau[ij] * r[i];
389  }
390 
391  sigma[j] = m_r[j].norm2(); // sigma[j] = r[j] * r[j];
392 
393  dcomplex r_0j = dotc(m_r[0], m_r[j]);
394  gamma_prime[j] = conj(r_0j) / sigma[j]; // gamma_prime[j] = (r[0] * r[j]) / sigma[j];
395 
396  //- a prescription to improve stability of BiCGStab(L)
397  double abs_rho = abs(r_0j) / sqrt(sigma[j] * sigma_0);
398  if (abs_rho < m_Omega_tolerance) {
399  gamma_prime[j] *= m_Omega_tolerance / abs_rho;
400  }
401  }
402 
403 
404  std::vector<dcomplex> gamma(m_N_L + 1);
405  gamma[N_L_tmp] = gamma_prime[N_L_tmp];
406 
407  for (int j = N_L_tmp - 1; j > 0; --j) {
408  dcomplex c_tmp = cmplx(0.0, 0.0);
409 
410  for (int i = j + 1; i < N_L_tmp + 1; ++i) {
411  int ji = index_ij(j, i);
412  c_tmp += tau[ji] * gamma[i];
413  }
414 
415  gamma[j] = gamma_prime[j] - c_tmp;
416  }
417 
418 
419  // NB. gamma_double_prime(m_N_L), not (m_N_L+1)
420  std::vector<dcomplex> gamma_double_prime(m_N_L);
421 
422  for (int j = 1; j < N_L_tmp; ++j) {
423  dcomplex c_tmp = cmplx(0.0, 0.0);
424 
425  for (int i = j + 1; i < N_L_tmp; ++i) {
426  int ji = index_ij(j, i);
427  c_tmp += tau[ji] * gamma[i + 1];
428  }
429 
430  gamma_double_prime[j] = gamma[j + 1] + c_tmp;
431  }
432 
433 
434  axpy(m_x, gamma[1], m_r[0]); // x += gamma[ 1] * r[ 0];
435  axpy(m_r[0], -gamma_prime[N_L_tmp], m_r[N_L_tmp]); // r[0] -= gamma_prime[N_L_tmp] * r[N_L_tmp];
436  axpy(m_u[0], -gamma[N_L_tmp], m_u[N_L_tmp]); // u[0] -= gamma[ N_L_tmp] * u[N_L_tmp];
437 
438  for (int j = 1; j < N_L_tmp; ++j) {
439  axpy(m_x, gamma_double_prime[j], m_r[j]); // x += gamma_double_prime[j] * r[j];
440  axpy(m_r[0], -gamma_prime[j], m_r[j]); // r[0] -= gamma_prime[ j] * r[j];
441  axpy(m_u[0], -gamma[j], m_u[j]); // u[0] -= gamma[ j] * u[j];
442  }
443 
444  rr = m_r[0].norm2(); // rr = r[0] * r[0];
445 
446 
447 #pragma omp barrier
448 #pragma omp master
449  {
450  m_rho_prev = rho_prev2;
451  m_alpha_prev = alpha_prev2;
452  m_N_L_prev = N_L_tmp;
453 
454  m_rho_prev *= -gamma_prime[N_L_tmp];
455  }
456 #pragma omp barrier
457 }
458 
459 
460 //====================================================================
462 {
463  const int NPE = CommonParameters::NPE();
464 
465  //- NB1 Nin = 2 * Nc * Nd, Nex = 1 for field_F
466  //- NB2 Nvol = CommonParameters::Nvol()/2 for eo
467  const int Nin = m_x.nin();
468  const int Nvol = m_x.nvol();
469  const int Nex = m_x.nex();
470 
471  const double gflop_fopr = m_fopr->flop_count();
472 
473  if (gflop_fopr < CommonParameters::epsilon_criterion()) {
474  vout.crucial(m_vl, "Warning at %s: no fopr->flop_count() is available, setting flop = 0\n", class_name.c_str());
475  return 0.0;
476  }
477 
478  const double gflop_axpy = (Nin * Nex * 2) * ((Nvol * NPE) / 1.0e+9);
479  const double gflop_dotc = (Nin * Nex * 4) * ((Nvol * NPE) / 1.0e+9);
480  const double gflop_norm = (Nin * Nex * 2) * ((Nvol * NPE) / 1.0e+9);
481 
482  const int N_L_prev_total = (m_Nconv_count - 1) / 2;
483 
484  const double gflop_init = gflop_fopr + gflop_axpy + gflop_norm;
485  const double gflop_step_BiCG_part = 2 * N_L_prev_total * gflop_fopr
486  + 3 * N_L_prev_total * gflop_dotc
487  + (N_L_prev_total + 2 * m_N_L_part_count) * gflop_axpy
488  + N_L_prev_total * gflop_norm;
489  const double gflop_step_L_part = (m_N_L_part_count + N_L_prev_total) * gflop_dotc
490  + (m_N_L_part_count + 3 * N_L_prev_total) * gflop_axpy
491  + (N_L_prev_total + m_Niter_count) * gflop_norm;
492  const double gflop_step = gflop_step_BiCG_part + gflop_step_L_part;
493  const double gflop_true_residual = gflop_fopr + gflop_axpy + gflop_norm;
494 
495  const double gflop = gflop_norm + gflop_init + gflop_step + gflop_true_residual * (m_Nrestart_count + 1)
496  + gflop_init * m_Nrestart_count;
497 
498 
499  return gflop;
500 }
501 
502 
503 //====================================================================
504 //============================================================END=====
void set_parameters_DS_L(const int N_L, const double Tol_L)
BridgeIO vout
Definition: bridgeIO.cpp:503
int fetch_bool(const string &key, bool &value) const
Definition: parameters.cpp:391
void detailed(const char *format,...)
Definition: bridgeIO.cpp:216
static double epsilon_criterion()
void solve_init(const Field &, double &)
double norm2() const
Definition: field.cpp:592
void general(const char *format,...)
Definition: bridgeIO.cpp:197
Container of Field-type object.
Definition: field.h:45
int fetch_double(const string &key, double &value) const
Definition: parameters.cpp:327
int nvol() const
Definition: field.h:127
void set_parameters_BiCGStab_series(const double Omega_tolerance)
Class for parameters.
Definition: parameters.h:46
void copy(Field &y, const Field &x)
copy(y, x): y = x
Definition: field.cpp:532
void solve(Field &solution, const Field &source, int &Nconv, double &diff)
int square_non_zero(const double v)
int nin() const
Definition: field.h:126
dcomplex dotc(const Field &y, const Field &x)
Definition: field.cpp:155
int fetch_int(const string &key, int &value) const
Definition: parameters.cpp:346
virtual double flop_count()
returns the flop in giga unit
Definition: fopr.h:120
void aypx(const double a, Field &y, const Field &x)
aypx(y, a, x): y := a * y + x
Definition: field.cpp:612
int nex() const
Definition: field.h:128
void paranoiac(const char *format,...)
Definition: bridgeIO.cpp:235
void reset(const int Nin, const int Nvol, const int Nex, const element_type cmpl=Element_type::COMPLEX)
Definition: field.h:95
void axpy(Field &y, const double a, const Field &x)
axpy(y, a, x): y := a * x + y
Definition: field.cpp:319
void crucial(const char *format,...)
Definition: bridgeIO.cpp:178
static const std::string class_name
virtual void mult(Field &, const Field &)=0
multiplies fermion operator to a given field (2nd argument)
int non_negative(const int v)
int index_ij(const int i, const int j)
void set_parameters(const Parameters &params)
static void assert_single_thread(const std::string &class_name)
assert currently running on single thread.
string get_string(const string &key) const
Definition: parameters.cpp:221
Bridge::VerboseLevel m_vl
Definition: solver.h:63
static VerboseLevel set_verbose_level(const std::string &str)
Definition: bridgeIO.cpp:131
int size() const
Definition: field.h:132