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