bicg.hpp 3.32 KB
Newer Older
1
// Software License for MTL
2
//
3
4
5
6
7
// Copyright (c) 2007 The Trustees of Indiana University.
//               2008 Dresden University of Technology and the Trustees of Indiana University.
//               2010 SimuNova UG (haftungsbeschränkt), www.simunova.com.
// All rights reserved.
// Authors: Peter Gottschling and Andrew Lumsdaine
8
//
9
// This file is part of the Matrix Template Library
10
//
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
// See also license.mtl.txt in the distribution.

#ifndef ITL_BICG_INCLUDE
#define ITL_BICG_INCLUDE

#include <complex>
#include <boost/numeric/itl/itl_fwd.hpp>
#include <boost/numeric/mtl/concept/collection.hpp>
#include <boost/numeric/mtl/operation/conj.hpp>
#include <boost/numeric/mtl/operation/resource.hpp>
#include <boost/numeric/mtl/interface/vpt.hpp>

namespace itl {

/// Bi-Conjugate Gradient
26
template < typename LinearOperator, typename Vector,
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
	   typename Preconditioner, typename Iteration >
int bicg(const LinearOperator &A, Vector &x, const Vector &b,
	 const Preconditioner &M, Iteration& iter)
{
    mtl::vampir_trace<7003> tracer;
    using mtl::conj;
    typedef typename mtl::Collection<Vector>::value_type Scalar;
    Scalar     rho_1(0), rho_2(0), alpha(0), beta(0);
    Vector     r(b - A * x), z(resource(x)), p(resource(x)), q(resource(x)),
 	       r_tilde(r), z_tilde(resource(x)), p_tilde(resource(x)), q_tilde(resource(x));

    while ( ! iter.finished(r)) {
	++iter;
	z= solve(M, r);
	z_tilde= adjoint_solve(M, r_tilde);
	rho_1= dot(z_tilde, z);

	if (rho_1 == 0.) return iter.fail(2, "bicg breakdown");
	if (iter.first()) {
	    p= z;
	    p_tilde= z_tilde;
	} else {
49
	    beta= rho_1 / rho_2;
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
	    p= z + beta * p;
	    p_tilde= z_tilde + conj(beta) * p_tilde;
	}

	q= A * p;
	q_tilde= adjoint(A) * p_tilde;
	alpha= rho_1 / dot(p_tilde, q);

	x+= alpha * p;
	r-= alpha * q;
	r_tilde-= conj(alpha) * q_tilde;

	rho_2= rho_1;
    }
    return iter;
}

/// Solver class for BiCG method; right preconditioner ignored (prints warning if not identity)
68
template < typename LinearOperator, typename Preconditioner= pc::identity<LinearOperator>,
69
70
71
72
73
	   typename RightPreconditioner= pc::identity<LinearOperator> >
class bicg_solver
{
  public:
    /// Construct solver from a linear operator; generate (left) preconditioner from it
74
    explicit bicg_solver(const LinearOperator& A) : A(A), L(A)
75
76
77
78
79
80
    {
	if (!pc::static_is_identity<RightPreconditioner>::value)
	    std::cerr << "Right Preconditioner ignored!" << std::endl;
    }

    /// Construct solver from a linear operator and (left) preconditioner
81
    bicg_solver(const LinearOperator& A, const Preconditioner& L) : A(A), L(L)
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
    {
	if (!pc::static_is_identity<RightPreconditioner>::value)
	    std::cerr << "Right Preconditioner ignored!" << std::endl;
    }

    /// Solve linear system approximately as specified by \p iter
    template < typename HilbertSpaceB, typename HilbertSpaceX, typename Iteration >
    int solve(const HilbertSpaceB& b, HilbertSpaceX& x, Iteration& iter) const
    {
	return bicg(A, x, b, L, iter);
    }

    /// Perform one bicg iteration on linear system
    template < typename HilbertSpaceB, typename HilbertSpaceX >
    int solve(const HilbertSpaceB& b, HilbertSpaceX& x) const
    {
	itl::basic_iteration<double> iter(x, 1, 0, 0);
	return bicg(A, x, b, L, iter);
    }
101

102
103
104
105
106
107
108
109
110
111
112
113
114
115
  private:
    const LinearOperator& A;
    Preconditioner        L;
};

} // namespace itl

#endif // ITL_BICG_INCLUDE