ug4
bicgstab.h
Go to the documentation of this file.
1 /*
2  * Copyright (c) 2010-2015: G-CSC, Goethe University Frankfurt
3  * Author: Andreas Vogel
4  *
5  * This file is part of UG4.
6  *
7  * UG4 is free software: you can redistribute it and/or modify it under the
8  * terms of the GNU Lesser General Public License version 3 (as published by the
9  * Free Software Foundation) with the following additional attribution
10  * requirements (according to LGPL/GPL v3 §7):
11  *
12  * (1) The following notice must be displayed in the Appropriate Legal Notices
13  * of covered and combined works: "Based on UG4 (www.ug4.org/license)".
14  *
15  * (2) The following notice must be displayed at a prominent place in the
16  * terminal output of covered works: "Based on UG4 (www.ug4.org/license)".
17  *
18  * (3) The following bibliography is recommended for citation and must be
19  * preserved in all covered files:
20  * "Reiter, S., Vogel, A., Heppner, I., Rupp, M., and Wittum, G. A massively
21  * parallel geometric multigrid solver on hierarchically distributed grids.
22  * Computing and visualization in science 16, 4 (2013), 151-164"
23  * "Vogel, A., Reiter, S., Rupp, M., Nägel, A., and Wittum, G. UG4 -- a novel
24  * flexible software system for simulating pde based models on high performance
25  * computers. Computing and visualization in science 16, 4 (2013), 165-179"
26  *
27  * This program is distributed in the hope that it will be useful,
28  * but WITHOUT ANY WARRANTY; without even the implied warranty of
29  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
30  * GNU Lesser General Public License for more details.
31  */
32 
33 #ifndef __H__UG__LIB_DISC__OPERATOR__LINEAR_OPERATOR__BICGSTAB__
34 #define __H__UG__LIB_DISC__OPERATOR__LINEAR_OPERATOR__BICGSTAB__
35 
36 #include <iostream>
37 #include <string>
38 #include <sstream>
39 
44 #ifdef UG_PARALLEL
46 #endif
47 
48 namespace ug{
49 
51 
65 template <typename TVector>
66 class BiCGStab
67  : public IPreconditionedLinearOperatorInverse<TVector>
68 {
69  public:
71  typedef TVector vector_type;
72 
75 
76  protected:
81 
82  public:
85  m_numRestarts(0), m_minOrtho(0.0)
86  {};
87 
89  : base_type ( spPrecond ),
90  m_numRestarts(0), m_minOrtho(0.0)
91  {}
92 
95  : base_type(spPrecond, spConvCheck),
96  m_numRestarts(0), m_minOrtho(0.0)
97  {};
98 
100  virtual const char* name() const {return "BiCGStab";}
101 
103  virtual bool supports_parallel() const
104  {
105  if(preconditioner().valid())
106  return preconditioner()->supports_parallel();
107  return true;
108  }
109 
110  // Solve J(u)*x = b, such that x = J(u)^{-1} b
112  {
113  LS_PROFILE_BEGIN(LS_ApplyReturnDefect);
114 
115  // check correct storage type in parallel
116  #ifdef UG_PARALLEL
117  if(!b.has_storage_type(PST_ADDITIVE) || !x.has_storage_type(PST_CONSISTENT)){
118  UG_LOG("BICGStab:b_storage_type"<<b.has_storage_type(PST_ADDITIVE)<<"\n");
119  UG_LOG("BICGStab:x_storage_type"<<x.has_storage_type(PST_CONSISTENT)<<"\n");
120  UG_THROW("BiCGStab: Inadequate storage format of Vectors.");
121  }
122  #endif
123 
124  // build defect: r := b - A*x
125  linear_operator()->apply_sub(b, x);
126  vector_type& r = b;
127 
128  // create vectors
129  SmartPtr<vector_type> spR = r.clone_without_values(); vector_type& r0 = *spR;
130  SmartPtr<vector_type> spP = r.clone_without_values(); vector_type& p = *spP;
131  SmartPtr<vector_type> spV = r.clone_without_values(); vector_type& v = *spV;
132  SmartPtr<vector_type> spT = r.clone_without_values(); vector_type& t = *spT;
133  SmartPtr<vector_type> spS = r.clone_without_values(); vector_type& s = *spS;
134  SmartPtr<vector_type> spQ = x.clone_without_values(); vector_type& q = *spQ;
135 
136  // prepare convergence check
138 
139  // compute start defect norm
140  convergence_check()->start(r);
141 
142  // convert b to unique (should already be unique due to norm calculation)
143  #ifdef UG_PARALLEL
144  if(!r.change_storage_type(PST_UNIQUE))
145  UG_THROW("BiCGStab: Cannot convert b to a vector with the 'unique' parallel storage type.");
146  #endif
147 
148  // needed variables
149  number rho = 1, alpha = 1, omega = 1, norm_r0 = 0.0;
150 
151  // restart flag (set to true at first run)
152  bool bRestart = true;
153 
154  write_debugXR(x, r, convergence_check()->step(), 'i');
155 
156  // Iteration loop
157  while(!convergence_check()->iteration_ended())
158  {
159  // check for restart based on fixed step number restart
160  if(m_numRestarts > 0 &&
161  (convergence_check()->step() % m_numRestarts == 0))
162  {
163  std::stringstream ss; ss <<
164  "Restarting: at every "<<m_numRestarts<<" Iterations";
165  convergence_check()->print_line(ss.str());
166  bRestart = true;
167  }
168 
169  // check if start values have to be set
170  if(bRestart)
171  {
172  // reset arbitrary vectors
173  r0 = r;
174 
175  // make r additive unique
176  #ifdef UG_PARALLEL
177  if(!r0.change_storage_type(PST_UNIQUE))
178  UG_THROW("BiCGStab: Cannot convert r to unique vector.");
179  #endif
180 
181  // set start vectors and alpha, omega:
182  // This will lead to p = r, since beta = 0.0
183  p = 0.0; alpha = 0.0;
184  v = 0.0; omega = 1.0;
185 
186  // set rhoOld to 1 (rhoOld = rho, see below)
187  rho = 1.0;
188 
189  // remember start norm
190  norm_r0 = convergence_check()->defect();
191 
192  // remove restart flag
193  bRestart = false;
194  }
195 
196  // remember current rho
197  const number rhoOld = rho;
198 
199  // Compute rho new
200  if (!r.size())
201  rho = 1.0;
202  else
203  rho = VecProd(r0, r);
204 
205  // check for restart compare (r, r0) > m_minOrtho * ||r|| ||r0||
206  const number norm_r = convergence_check()->defect();
207  if(fabs(rho)/(norm_r * norm_r0) <= m_minOrtho)
208  {
209  std::stringstream ss; ss <<
210  "Restarting: Min Orthogonality "<<m_minOrtho<<" missed: "
211  <<"(r,r0)="<<fabs(rho)<<", ||r||="<<norm_r<<", ||r0||= "
212  <<norm_r0;
213  convergence_check()->print_line(ss.str());
214  bRestart = true;
215  }
216 
217  // check that rhoOld valid
218  if(rhoOld == 0.0)
219  {
220  UG_LOG("BiCGStab: Method breakdown with rhoOld = "<<rhoOld<<
221  ". Aborting iteration.\n");
222  return false;
223  }
224 
225  // Compute new beta
226  const number beta = (rho/rhoOld) * (alpha/omega);
227 
228  // update p = r + beta * p - beta * omega * v
229  VecScaleAdd(p, 1.0, r, beta, p, -beta*omega, v);
230 
231  // apply q = M^-1 * p ...
232  if(preconditioner().valid())
233  {
235  if(!preconditioner()->apply(q, p))
236  {
237  UG_LOG("BiCGStab: Cannot apply preconditioner. Aborting.\n");
239  return false;
240  }
242  }
243  // ... or copy q = p
244  else
245  {
246  q = p;
247 
248  // make q consistent
249  #ifdef UG_PARALLEL
250  if(!q.change_storage_type(PST_CONSISTENT))
251  UG_THROW("BiCGStab: Cannot convert q to consistent vector.");
252  #endif
253  }
254 
255  // post-process the correction
257 
258  // compute v := A*q
259  linear_operator()->apply(v, q);
260 
261  // make v unique
262  #ifdef UG_PARALLEL
263  if(!v.change_storage_type(PST_UNIQUE))
264  UG_THROW("BiCGStab: Cannot convert v to unique vector.");
265  #endif
266 
267  // alpha = (v,r)
268  if (!v.size())
269  alpha = 1.0;
270  else
271  alpha = VecProd(v, r0);
272 
273  // check validity of alpha
274  if(alpha == 0.0){
275  UG_LOG("BiCGStab: Method breakdown: alpha = "<<alpha<<
276  " is an invalid value. Aborting iteration.\n");
277  return false;
278  }
279 
280  // alpha = rho/(v,r)
281  alpha = rho/alpha;
282 
283  // add: x := x + alpha * q
284  VecScaleAdd(x, 1.0, x, alpha, q);
285 
286  // compute s = r - alpha*v
287  VecScaleAdd(s, 1.0, r, -alpha, v);
288 
289  // check convergence
290  convergence_check()->update(s);
291 
292  write_debugXR(x, s, convergence_check()->step(), 'a');
293 
294  // if finished: set output to last defect and exist loop
295  if(convergence_check()->iteration_ended())
296  {
297  r = s; break;
298  }
299 
300  // apply q = M^-1 * s ...
301  if(preconditioner().valid())
302  {
304  if(!preconditioner()->apply(q, s))
305  {
306  UG_LOG("BiCGStab: Cannot apply preconditioner. Aborting.\n");
308  return false;
309  }
311  }
312  // ... or set q:=s
313  else
314  {
315  q = s;
316 
317  // make q consistent
318  #ifdef UG_PARALLEL
319  if(!q.change_storage_type(PST_CONSISTENT))
320  UG_THROW("BiCGStab: Cannot convert q to consistent vector.");
321  #endif
322  }
323 
324  // post-process the correction
326 
327  // compute t := A*q
328  linear_operator()->apply(t, q);
329 
330  // make t unique
331  #ifdef UG_PARALLEL
332  if(!t.change_storage_type(PST_UNIQUE))
333  UG_THROW("BiCGStab: Cannot convert t to unique vector.");
334  #endif
335 
336  // tt = (t,t)
337  number tt;
338  if (!t.size())
339  tt = 1.0;
340  else
341  tt = VecProd(t, t);
342 
343  // omega = (s,t)
344  if (!s.size())
345  omega = 1.0;
346  else
347  omega = VecProd(s, t);
348 
349  // check tt
350  if(tt == 0.0)
351  {
352  UG_LOG("BiCGStab: Method breakdown tt = "<<tt<<" is an "
353  "invalid value. Aborting iteration.\n");
354  return false;
355  }
356 
357  // omega = (s,t)/(t,t)
358  omega = omega/tt;
359 
360  // add: x := x + omega * q
361  VecScaleAdd(x, 1.0, x, omega, q);
362 
363  // compute r = s - omega*t
364  VecScaleAdd(r, 1.0, s, -omega, t);
365 
366  // check convergence
367  convergence_check()->update(r);
368 
369  write_debugXR(x, r, convergence_check()->step(), 'b');
370 
371  // check values
372  if(omega == 0.0)
373  {
374  UG_LOG("BiCGStab: Method breakdown with omega = "<<omega<<
375  ". Aborting iteration.\n");
376  return false;
377  }
378  }
379 
380  // print ending output
381  return convergence_check()->post();
382  }
383 
384 
386  void set_restart(int numRestarts) {m_numRestarts = numRestarts;}
387 
389  void set_min_orthogonality(number minOrtho) {m_minOrtho = minOrtho;}
390 
393  {
395  }
396 
399  {
401  }
402 
403  protected:
406  {
407  // set iteration symbol and name
408  convergence_check()->set_name(name());
409  convergence_check()->set_symbol('%');
410 
411  // set preconditioner string
412  std::string s;
413  if(preconditioner().valid())
414  s = std::string(" (Precond: ") + preconditioner()->name() + ")";
415  else
416  s = " (No Preconditioner) ";
417  convergence_check()->set_info(s);
418  }
419 
421  void write_debugXR(vector_type &x, vector_type &r, int loopCnt, char phase)
422  {
423  if(!this->vector_debug_writer_valid()) return;
424  char ext[20]; sprintf(ext, "-%c_iter%03d", phase, loopCnt);
425  write_debug(r, std::string("BiCGStab_Residual") + ext + ".vec");
426  write_debug(x, std::string("BiCGStab_Solution") + ext + ".vec");
427  }
428 
430  void enter_precond_debug_section(int loopCnt, char phase)
431  {
432  if(!this->vector_debug_writer_valid()) return;
433  char ext[20]; sprintf(ext, "-%c_iter%03d", phase, loopCnt);
434  this->enter_vector_debug_writer_section(std::string("BiCGStab_Precond") + ext);
435  }
436 
437  public:
438  virtual std::string config_string() const
439  {
440  std::stringstream ss;
441  ss << "BiCGStab( restart = " << m_numRestarts << ", min_orthogonality = " << m_minOrtho << ")\n";
443  return ss.str();
444  }
445  protected:
448 
451 
453 
459 };
460 
461 } // end namespace ug
462 
463 #endif /* __H__UG__LIB_DISC__OPERATOR__LINEAR_OPERATOR__BICGSTAB__ */
parameterString p
parameterString s
Definition: smart_pointer.h:108
the BiCGStab method as a solver for linear operators
Definition: bicgstab.h:68
void add_postprocess_corr(SmartPtr< IPProcessVector< vector_type > > p)
adds a post-process for the iterates
Definition: bicgstab.h:392
BiCGStab()
constructors
Definition: bicgstab.h:84
void write_debugXR(vector_type &x, vector_type &r, int loopCnt, char phase)
debugger output: solution and residual
Definition: bicgstab.h:421
void enter_precond_debug_section(int loopCnt, char phase)
debugger section for the preconditioner
Definition: bicgstab.h:430
BiCGStab(SmartPtr< ILinearIterator< vector_type, vector_type > > spPrecond)
Definition: bicgstab.h:88
TVector vector_type
Vector type.
Definition: bicgstab.h:71
void remove_postprocess_corr(SmartPtr< IPProcessVector< vector_type > > p)
removes a post-process for the iterates
Definition: bicgstab.h:398
void set_restart(int numRestarts)
sets to restart at given number of iteration steps
Definition: bicgstab.h:386
void set_min_orthogonality(number minOrtho)
sets to restart if given orthogonality missed
Definition: bicgstab.h:389
BiCGStab(SmartPtr< ILinearIterator< vector_type > > spPrecond, SmartPtr< IConvergenceCheck< vector_type > > spConvCheck)
Definition: bicgstab.h:93
virtual const char * name() const
name of solver
Definition: bicgstab.h:100
number m_minOrtho
minimal value in (0,1) accepted for Orthoginality before restart
Definition: bicgstab.h:450
SmartPtr< ILinearIterator< X, X > > preconditioner()
Definition: preconditioned_linear_operator_inverse.h:118
PProcessChain< vector_type > m_corr_post_process
postprocessor for the correction in the iterations
Definition: bicgstab.h:458
IPreconditionedLinearOperatorInverse< vector_type > base_type
Base type.
Definition: bicgstab.h:74
void prepare_conv_check()
prepares the output of the convergence check
Definition: bicgstab.h:405
virtual bool supports_parallel() const
returns if parallel solving is supported
Definition: bicgstab.h:103
virtual std::string config_string() const
returns information about configuration parameters
Definition: bicgstab.h:438
virtual bool apply_return_defect(vector_type &x, vector_type &b)
Definition: bicgstab.h:111
int m_numRestarts
restarts at every numRestarts steps (numRestarts <= 0 --> never)
Definition: bicgstab.h:447
Definition: convergence_check.h:72
describes a linear iterator
Definition: linear_iterator.h:81
ConstSmartPtr< IConvergenceCheck< X > > convergence_check() const
returns the convergence check
Definition: linear_operator_inverse.h:202
SmartPtr< ILinearOperator< Y, X > > linear_operator()
returns the current Operator this Inverse Operator is initialized for
Definition: linear_operator_inverse.h:227
int step() const
returns the current number of steps
Definition: linear_operator_inverse.h:211
interface for pre- and postprocess functions
Definition: pprocess.h:59
describes an inverse linear mapping X->X
Definition: preconditioned_linear_operator_inverse.h:63
virtual bool apply(TVector &x, const TVector &b)
Definition: preconditioned_linear_operator_inverse.h:152
SmartPtr< ILinearIterator< X, X > > preconditioner()
Definition: preconditioned_linear_operator_inverse.h:118
std::string config_string_preconditioner_convergence_check() const
returns config information of convergence check and preconditioner
Definition: preconditioned_linear_operator_inverse.h:200
void add(SmartPtr< p_process_type > p)
adds an operation at the end of the chain
Definition: pprocess.h:96
void apply(vector_type &v)
performs all the operations
Definition: pprocess.h:86
void remove(SmartPtr< p_process_type > p)
removes an operation from the chain
Definition: pprocess.h:108
bool vector_debug_writer_valid() const
returns true if the debug writer is set
Definition: debug_writer.h:290
void enter_vector_debug_writer_section(std::string secDir)
enters a debugging section
Definition: debug_writer.h:331
void leave_vector_debug_writer_section()
leaves a debugging section
Definition: debug_writer.h:342
void write_debug(const vector_type &vec, const char *filename)
writing debug output for a vector (if debug writer set)
Definition: debug_writer.h:293
TVector vector_type
type of vector
Definition: debug_writer.h:269
@ PST_CONSISTENT
Definition: parallel_storage_type.h:68
@ PST_UNIQUE
Definition: parallel_storage_type.h:70
@ PST_ADDITIVE
Definition: parallel_storage_type.h:69
function common sprintf(s,...)
number omega
number alpha
#define UG_THROW(msg)
Definition: error.h:57
#define UG_LOG(msg)
Definition: log.h:367
double number
Definition: types.h:124
void VecScaleAdd(vector_t &vOut, typename vector_t::value_type s1, const vector_t &v1, typename vector_t::value_type s2, const vector_t &v2)
Scales two Vectors, adds them and returns the sum in a third vector.
Definition: math_vector_functions_common_impl.hpp:265
#define LS_PROFILE_BEGIN(name)
Definition: linear_solver_profiling.h:40
the ug namespace
double VecProd(const double &a, const double &b)
returns scal<a, b>
Definition: operations_vec.h:84