Loading [MathJax]/extensions/tex2jax.js
ug4
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
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
48
49namespace ug{
50
52
66template <typename TVector>
69{
70 public:
72 typedef TVector vector_type;
73
76
77 protected:
82
83 public:
87 {};
88
93
96 : base_type(spPrecond, spConvCheck),
98 {};
99
101 virtual const char* name() const {return "BiCGStab";}
102
104 virtual bool supports_parallel() const
105 {
106 if(preconditioner().valid())
107 return preconditioner()->supports_parallel();
108 return true;
109 }
110
111 // Solve J(u)*x = b, such that x = J(u)^{-1} b
113 {
114 LS_PROFILE_BEGIN(LS_ApplyReturnDefect);
115
116 // check correct storage type in parallel
117 #ifdef UG_PARALLEL
118 if(!b.has_storage_type(PST_ADDITIVE) || !x.has_storage_type(PST_CONSISTENT)){
119 UG_LOG("BICGStab:b_storage_type"<<b.has_storage_type(PST_ADDITIVE)<<"\n");
120 UG_LOG("BICGStab:x_storage_type"<<x.has_storage_type(PST_CONSISTENT)<<"\n");
121 UG_THROW("BiCGStab: Inadequate storage format of Vectors.");
122 }
123 #endif
124
125 // build defect: r := b - A*x
126 linear_operator()->apply_sub(b, x);
127 vector_type& r = b;
128
129 // create vectors
130 SmartPtr<vector_type> spR = r.clone_without_values(); vector_type& r0 = *spR;
131 SmartPtr<vector_type> spP = r.clone_without_values(); vector_type& p = *spP;
132 SmartPtr<vector_type> spV = r.clone_without_values(); vector_type& v = *spV;
133 SmartPtr<vector_type> spT = r.clone_without_values(); vector_type& t = *spT;
134 SmartPtr<vector_type> spS = r.clone_without_values(); vector_type& s = *spS;
135 SmartPtr<vector_type> spQ = x.clone_without_values(); vector_type& q = *spQ;
136
137 // prepare convergence check
139
140 // compute start defect norm
141 convergence_check()->start(r);
142
143 // convert b to unique (should already be unique due to norm calculation)
144 #ifdef UG_PARALLEL
145 if(!r.change_storage_type(PST_UNIQUE))
146 UG_THROW("BiCGStab: Cannot convert b to a vector with the 'unique' parallel storage type.");
147 #endif
148
149 // needed variables
150 number rho = 1, alpha = 1, omega = 1, norm_r0 = 0.0;
151
152 // restart flag (set to true at first run)
153 bool bRestart = true;
154
155 write_debugXR(x, r, convergence_check()->step(), 'i');
156
157 // Iteration loop
158 while(!convergence_check()->iteration_ended())
159 {
160 // check for restart based on fixed step number restart
161 if(m_numRestarts > 0 &&
163 {
164 std::stringstream ss; ss <<
165 "Restarting: at every "<<m_numRestarts<<" Iterations";
166 convergence_check()->print_line(ss.str());
167 bRestart = true;
168 }
169
170 // check if start values have to be set
171 if(bRestart)
172 {
173 // reset arbitrary vectors
174 r0 = r;
175
176 // make r additive unique
177 #ifdef UG_PARALLEL
178 if(!r0.change_storage_type(PST_UNIQUE))
179 UG_THROW("BiCGStab: Cannot convert r to unique vector.");
180 #endif
181
182 // set start vectors and alpha, omega:
183 // This will lead to p = r, since beta = 0.0
184 p = 0.0; alpha = 0.0;
185 v = 0.0; omega = 1.0;
186
187 // set rhoOld to 1 (rhoOld = rho, see below)
188 rho = 1.0;
189
190 // remember start norm
191 norm_r0 = convergence_check()->defect();
192
193 // remove restart flag
194 bRestart = false;
195 }
196
197 // remember current rho
198 const number rhoOld = rho;
199
200 // Compute rho new
201 if (!r.size())
202 rho = 1.0;
203 else
204 rho = VecProd(r0, r);
205
206 // check for restart compare (r, r0) > m_minOrtho * ||r|| ||r0||
207 const number norm_r = convergence_check()->defect();
208 if(fabs(rho)/(norm_r * norm_r0) <= m_minOrtho)
209 {
210 std::stringstream ss; ss <<
211 "Restarting: Min Orthogonality "<<m_minOrtho<<" missed: "
212 <<"(r,r0)="<<fabs(rho)<<", ||r||="<<norm_r<<", ||r0||= "
213 <<norm_r0;
214 convergence_check()->print_line(ss.str());
215 bRestart = true;
216 }
217
218 // check that rhoOld valid
219 if(rhoOld == 0.0)
220 {
221 UG_LOG("BiCGStab: Method breakdown with rhoOld = "<<rhoOld<<
222 ". Aborting iteration.\n");
223 return false;
224 }
225
226 // Compute new beta
227 const number beta = (rho/rhoOld) * (alpha/omega);
228
229 // update p = r + beta * p - beta * omega * v
230 VecScaleAdd(p, 1.0, r, beta, p, -beta*omega, v);
231
232 // apply q = M^-1 * p ...
233 if(preconditioner().valid())
234 {
236 if(!preconditioner()->apply(q, p))
237 {
238 UG_LOG("BiCGStab: Cannot apply preconditioner. Aborting.\n");
240 return false;
241 }
243 }
244 // ... or copy q = p
245 else
246 {
247 q = p;
248
249 // make q consistent
250 #ifdef UG_PARALLEL
251 if(!q.change_storage_type(PST_CONSISTENT))
252 UG_THROW("BiCGStab: Cannot convert q to consistent vector.");
253 #endif
254 }
255
256 // post-process the correction
258
259 // compute v := A*q
260 linear_operator()->apply(v, q);
261
262 // make v unique
263 #ifdef UG_PARALLEL
264 if(!v.change_storage_type(PST_UNIQUE))
265 UG_THROW("BiCGStab: Cannot convert v to unique vector.");
266 #endif
267
268 // alpha = (v,r)
269 if (!v.size())
270 alpha = 1.0;
271 else
272 alpha = VecProd(v, r0);
273
274 // check validity of alpha
275 if(alpha == 0.0){
276 UG_LOG("BiCGStab: Method breakdown: alpha = "<<alpha<<
277 " is an invalid value. Aborting iteration.\n");
278 return false;
279 }
280
281 // alpha = rho/(v,r)
282 alpha = rho/alpha;
283
284 // add: x := x + alpha * q
285 VecScaleAdd(x, 1.0, x, alpha, q);
286
287 // compute s = r - alpha*v
288 VecScaleAdd(s, 1.0, r, -alpha, v);
289
290 // check convergence
291 convergence_check()->update(s);
292
293 write_debugXR(x, s, convergence_check()->step(), 'a');
294
295 // if finished: set output to last defect and exist loop
296 if(convergence_check()->iteration_ended())
297 {
298 r = s; break;
299 }
300
301 // apply q = M^-1 * s ...
302 if(preconditioner().valid())
303 {
305 if(!preconditioner()->apply(q, s))
306 {
307 UG_LOG("BiCGStab: Cannot apply preconditioner. Aborting.\n");
309 return false;
310 }
312 }
313 // ... or set q:=s
314 else
315 {
316 q = s;
317
318 // make q consistent
319 #ifdef UG_PARALLEL
320 if(!q.change_storage_type(PST_CONSISTENT))
321 UG_THROW("BiCGStab: Cannot convert q to consistent vector.");
322 #endif
323 }
324
325 // post-process the correction
327
328 // compute t := A*q
329 linear_operator()->apply(t, q);
330
331 // make t unique
332 #ifdef UG_PARALLEL
333 if(!t.change_storage_type(PST_UNIQUE))
334 UG_THROW("BiCGStab: Cannot convert t to unique vector.");
335 #endif
336
337 // tt = (t,t)
338 number tt;
339 if (!t.size())
340 tt = 1.0;
341 else
342 tt = VecProd(t, t);
343
344 // omega = (s,t)
345 if (!s.size())
346 omega = 1.0;
347 else
348 omega = VecProd(s, t);
349
350 // check tt
351 if(tt == 0.0)
352 {
353 UG_LOG("BiCGStab: Method breakdown tt = "<<tt<<" is an "
354 "invalid value. Aborting iteration.\n");
355 return false;
356 }
357
358 // omega = (s,t)/(t,t)
359 omega = omega/tt;
360
361 // add: x := x + omega * q
362 VecScaleAdd(x, 1.0, x, omega, q);
363
364 // compute r = s - omega*t
365 VecScaleAdd(r, 1.0, s, -omega, t);
366
367 // check convergence
368 convergence_check()->update(r);
369
370 write_debugXR(x, r, convergence_check()->step(), 'b');
371
372 // check values
373 if(omega == 0.0)
374 {
375 UG_LOG("BiCGStab: Method breakdown with omega = "<<omega<<
376 ". Aborting iteration.\n");
377 return false;
378 }
379 }
380
381 // print ending output
382 return convergence_check()->post();
383 }
384
385
387 void set_restart(int numRestarts) {m_numRestarts = numRestarts;}
388
390 void set_min_orthogonality(number minOrtho) {m_minOrtho = minOrtho;}
391
397
403
404 protected:
407 {
408 // set iteration symbol and name
409 convergence_check()->set_name(name());
410 convergence_check()->set_symbol('%');
411
412 // set preconditioner string
413 std::string s;
414 if(preconditioner().valid())
415 s = std::string(" (Precond: ") + preconditioner()->name() + ")";
416 else
417 s = " (No Preconditioner) ";
418 convergence_check()->set_info(s);
419 }
420
422 void write_debugXR(vector_type &x, vector_type &r, int loopCnt, char phase)
423 {
424 if(!this->vector_debug_writer_valid()) return;
425 std::string ext = GetStringPrintf("-%c_iter%03d", phase, loopCnt);
426 write_debug(r, std::string("BiCGStab_Residual") + ext + ".vec");
427 write_debug(x, std::string("BiCGStab_Solution") + ext + ".vec");
428 }
429
431 void enter_precond_debug_section(int loopCnt, char phase)
432 {
433 if(!this->vector_debug_writer_valid()) return;
434 std::string ext = GetStringPrintf("-%c_iter%03d", phase, loopCnt);
435 this->enter_vector_debug_writer_section(std::string("BiCGStab_Precond") + ext);
436 }
437
438 public:
439 virtual std::string config_string() const
440 {
441 std::stringstream ss;
442 ss << "BiCGStab( restart = " << m_numRestarts << ", min_orthogonality = " << m_minOrtho << ")\n";
444 return ss.str();
445 }
446 protected:
449
452
454
460};
461
462} // end namespace ug
463
464#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:69
void add_postprocess_corr(SmartPtr< IPProcessVector< vector_type > > p)
adds a post-process for the iterates
Definition bicgstab.h:393
BiCGStab()
constructors
Definition bicgstab.h:85
void write_debugXR(vector_type &x, vector_type &r, int loopCnt, char phase)
debugger output: solution and residual
Definition bicgstab.h:422
void enter_precond_debug_section(int loopCnt, char phase)
debugger section for the preconditioner
Definition bicgstab.h:431
BiCGStab(SmartPtr< ILinearIterator< vector_type, vector_type > > spPrecond)
Definition bicgstab.h:89
TVector vector_type
Vector type.
Definition bicgstab.h:72
void remove_postprocess_corr(SmartPtr< IPProcessVector< vector_type > > p)
removes a post-process for the iterates
Definition bicgstab.h:399
void set_restart(int numRestarts)
sets to restart at given number of iteration steps
Definition bicgstab.h:387
void set_min_orthogonality(number minOrtho)
sets to restart if given orthogonality missed
Definition bicgstab.h:390
virtual const char * name() const
name of solver
Definition bicgstab.h:101
BiCGStab(SmartPtr< ILinearIterator< vector_type > > spPrecond, SmartPtr< IConvergenceCheck< vector_type > > spConvCheck)
Definition bicgstab.h:94
number m_minOrtho
minimal value in (0,1) accepted for Orthoginality before restart
Definition bicgstab.h:451
PProcessChain< vector_type > m_corr_post_process
postprocessor for the correction in the iterations
Definition bicgstab.h:459
IPreconditionedLinearOperatorInverse< vector_type > base_type
Base type.
Definition bicgstab.h:75
void prepare_conv_check()
prepares the output of the convergence check
Definition bicgstab.h:406
virtual bool supports_parallel() const
returns if parallel solving is supported
Definition bicgstab.h:104
virtual std::string config_string() const
returns information about configuration parameters
Definition bicgstab.h:439
virtual bool apply_return_defect(vector_type &x, vector_type &b)
Definition bicgstab.h:112
int m_numRestarts
restarts at every numRestarts steps (numRestarts <= 0 --> never)
Definition bicgstab.h:448
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
int step() const
returns the current number of steps
Definition linear_operator_inverse.h:211
SmartPtr< ILinearOperator< Y, X > > linear_operator()
returns the current Operator this Inverse Operator is initialized for
Definition linear_operator_inverse.h:227
interface for pre- and postprocess functions
Definition pprocess.h:59
describes an inverse linear mapping X->X
Definition preconditioned_linear_operator_inverse.h:63
SmartPtr< ILinearIterator< vector_type, vector_type > > preconditioner()
Definition preconditioned_linear_operator_inverse.h:118
virtual bool apply(TVector &x, const TVector &b)
Definition preconditioned_linear_operator_inverse.h:152
std::string config_string_preconditioner_convergence_check() const
returns config information of convergence check and preconditioner
Definition preconditioned_linear_operator_inverse.h:200
a chain of pre- or postprocess operations
Definition pprocess.h:78
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
@ PST_CONSISTENT
Definition parallel_storage_type.h:68
@ PST_UNIQUE
Definition parallel_storage_type.h:70
@ PST_ADDITIVE
Definition parallel_storage_type.h:69
#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
std::string GetStringPrintf(const char *format, Args... args)
Definition string_util.h:449