26   bool init = Solver::Factory::Register(
"BiCGStab_IDS_L_Cmplx", create_object);
 
   44 #ifdef USE_PARAMETERS_FACTORY 
   59   const string str_vlevel = params.
get_string(
"verbose_level");
 
   70   err += params.
fetch_int(
"maximum_number_of_iteration", Niter);
 
   71   err += params.
fetch_double(
"convergence_criterion_squared", Stop_cond);
 
   72   err += params.
fetch_int(
"number_of_orthonormal_vectors", N_L);
 
   73   err += params.
fetch_double(
"tolerance_for_DynamicSelection_of_L", Tol_L);
 
  137                                         int& Nconv, 
double& diff)
 
  139   double bnorm2 = b.
norm2();
 
  140   double snorm  = 1.0 / bnorm2;
 
  141   int    bsize  = b.
size();
 
  157   bool is_converged = 
false;
 
  162   for (
int iter = 0; iter < 
m_Niter; iter++) {
 
  163     if (is_converged) 
break;
 
  222     if ((
s.
nin() != Nin) || (
s.
nvol() != Nvol) || (
s.
nex() != Nex)) {
 
  232     for (
int i = 0; i < 
m_N_L + 1; ++i) {
 
  233       u[i].reset(Nin, Nvol, Nex);
 
  234       r[i].reset(Nin, Nvol, Nex);
 
  279   dcomplex c_Rayleigh_prev = cmplx(0.0, 0.0);
 
  281   bool is_converged_L = 
false;
 
  284   for (
int j = 0; j < 
N_L_prev; ++j) {
 
  285     if (!is_converged_L) {
 
  289       dcomplex beta = alpha_prev2 * (rho / rho_prev2);
 
  293       for (
int i = 0; i < j + 1; ++i) {
 
  301       alpha_prev2 = rho_prev2 / conj(gamma);
 
  303       for (
int i = 0; i < j + 1; ++i) {
 
  304         axpy(
r[i], -alpha_prev2, 
u[i + 1]);  
 
  309       axpy(
x, alpha_prev2, 
u[0]);          
 
  315       double   r_tmp      = 
r[j].norm2();                 
 
  316       dcomplex c_Rayleigh = 
dotc(
r[j], 
r[j + 1]) / r_tmp; 
 
  318       dcomplex c_E = (c_Rayleigh - c_Rayleigh_prev) / c_Rayleigh;
 
  321       c_Rayleigh_prev = c_Rayleigh;
 
  331         is_converged_L = 
true;
 
  337   std::vector<double>   sigma(
m_N_L + 1);
 
  338   std::vector<dcomplex> gamma_prime(
m_N_L + 1);
 
  341   std::vector<dcomplex> tau(
m_N_L * (
m_N_L + 1));
 
  344   for (
int j = 1; j < N_L_tmp + 1; ++j) {
 
  345     for (
int i = 1; i < j; ++i) {
 
  348       dcomplex r_ji = 
dotc(
r[j], 
r[i]);
 
  349       tau[ij] = conj(r_ji) / sigma[i];  
 
  350       axpy(
r[j], -tau[ij], 
r[i]);       
 
  353     sigma[j] = 
r[j].norm2();  
 
  355     dcomplex r_0j = 
dotc(
r[0], 
r[j]);
 
  356     gamma_prime[j] = conj(r_0j) / sigma[j]; 
 
  360   std::vector<dcomplex> gamma(
m_N_L + 1);
 
  363   gamma[N_L_tmp] = gamma_prime[N_L_tmp];
 
  365   for (
int j = N_L_tmp - 1; j > 0; --j) {
 
  366     c_tmp = cmplx(0.0, 0.0);
 
  368     for (
int i = j + 1; i < N_L_tmp + 1; ++i) {
 
  370       c_tmp += tau[ji] * gamma[i];
 
  373     gamma[j] = gamma_prime[j] - c_tmp;
 
  378   std::vector<dcomplex> gamma_double_prime(
m_N_L);
 
  380   for (
int j = 1; j < N_L_tmp; ++j) {
 
  381     c_tmp = cmplx(0.0, 0.0);
 
  383     for (
int i = j + 1; i < N_L_tmp; ++i) {
 
  385       c_tmp += tau[ji] * gamma[i + 1];
 
  388     gamma_double_prime[j] = gamma[j + 1] + c_tmp;
 
  393   axpy(
r[0], -gamma_prime[N_L_tmp], 
r[N_L_tmp]);  
 
  394   axpy(
u[0], -gamma[N_L_tmp], 
u[N_L_tmp]);        
 
  396   for (
int j = 1; j < N_L_tmp; ++j) {
 
  397     axpy(
x, gamma_double_prime[j], 
r[j]);      
 
  398     axpy(
r[0], -gamma_prime[j], 
r[j]);         
 
  399     axpy(
u[0], -gamma[j], 
u[j]);               
 
int index_ij(int i, int j)
 
void detailed(const char *format,...)
 
void Register_string(const string &, const string &)
 
void set_parameters_DS_L(const int N_L, const double Tol_L)
 
void general(const char *format,...)
 
void Register_int(const string &, const int)
 
Container of Field-type object. 
 
void copy(Field &y, const Field &x)
copy(y, x): y = x 
 
BiCGStab(IDS_L) algorithm. 
 
void solve_init(const Field &, double &)
 
Parameters_Solver_BiCGStab_IDS_L_Cmplx()
 
int square_non_zero(const double v)
 
dcomplex dotc(const Field &y, const Field &x)
 
void reset(const int Nin, const int Nvol, const int Nex, const element_type cmpl=COMPLEX)
 
void aypx(const double a, Field &y, const Field &x)
aypx(y, a, x): y := a * y + x 
 
void paranoiac(const char *format,...)
 
void axpy(Field &y, const double a, const Field &x)
axpy(y, a, x): y := a * x + y 
 
static const std::string class_name
 
void crucial(const char *format,...)
 
Base class for linear solver class family. 
 
static bool Register(const std::string &realm, const creator_callback &cb)
 
virtual void mult(Field &, const Field &)=0
multiplies fermion operator to a given field (2nd argument) 
 
void solve_step(double &)
 
void reset_field(const Field &)
 
int non_negative(const int v)
 
void Register_double(const string &, const double)
 
void set_parameters(const Parameters ¶ms)
 
Base class of fermion operator family. 
 
int fetch_double(const string &key, double &val) const 
 
string get_string(const string &key) const 
 
void solve(Field &solution, const Field &source, int &Nconv, double &diff)
 
int fetch_int(const string &key, int &val) const 
 
Bridge::VerboseLevel m_vl
 
static VerboseLevel set_verbose_level(const std::string &str)
 
static void assert_single_thread(const std::string &classname)
assert currently running on single thread.