Loading [MathJax]/extensions/tex2jax.js

DUNE-FEM (unstable)

cg.hh
1#ifndef DUNE_FEM_CG_HH
2#define DUNE_FEM_CG_HH
3
4#include <type_traits>
5#include <vector>
6#include <cmath>
7
9#include <dune/fem/solver/parameter.hh>
10
11namespace Dune
12{
13namespace Fem
14{
15namespace LinearSolver
16{
17
18 template <class Operator, class Precoditioner, class DiscreteFunction>
19 inline int cg( Operator &op,
20 Precoditioner* preconditioner,
21 std::vector< DiscreteFunction >& tempMem,
22 DiscreteFunction& x,
23 const DiscreteFunction& b,
24 const double epsilon,
25 const int maxIterations,
26 const int toleranceCriteria,
27 std::ostream* os = nullptr )
28 {
29 typedef typename DiscreteFunction::RangeFieldType RangeFieldType;
30 typedef typename Dune::FieldTraits< RangeFieldType >::real_type RealType;
31
32 assert( preconditioner ? tempMem.size() == 5 : tempMem.size() == 3 );
33
34 DiscreteFunction& h = tempMem[ 0 ];
35
36 //h=Ax
37 op( x, h );
38
39 //r=Ax-b
40 DiscreteFunction& r = tempMem[ 1 ];
41 r.assign( h );
42 r -= b;
43
44 //p=b-A*x <= r_0 Deufelhard
45 DiscreteFunction& p = tempMem[ 2 ];
46 p.assign( b );
47 p -= h;
48
49 DiscreteFunction& s = ( preconditioner ) ? tempMem[ 3 ] : p;
50
51 //q=B*p <=q Deufelhard
52 DiscreteFunction& q = ( preconditioner ) ? tempMem[ 4 ] : p;
53 if( preconditioner )
54 {
55 (*preconditioner)( p, q );
56 s.assign( q );
57 }
58
59 RangeFieldType prevResidual = 0; // note that these will be real_type but require scalar product evaluation
60 RangeFieldType residual = p.scalarProductDofs( q );//<p,Bp>
61
62 const RealType tolerance = epsilon * epsilon * (
63 toleranceCriteria == ToleranceCriteria::relative ? b.normSquaredDofs( ) :
64 toleranceCriteria == ToleranceCriteria::residualReduction ? std::real(residual) :
65 // absolute tolerance
66 1.0
67 );
68
69 int iterations = 0;
70 for( iterations = 0; (std::real(residual) > tolerance) && (iterations < maxIterations); ++iterations )
71 {
72 if( iterations > 0 )
73 {
74 assert( residual/prevResidual == residual/prevResidual );
75 const RangeFieldType beta= residual / prevResidual;
76 q *= beta;
77 if( preconditioner )
78 {
79 q += s;
80 }
81 else
82 {
83 p -= r;
84 }
85 }
86
87 op( q, h );
88
89 RangeFieldType qdoth = q.scalarProductDofs( h );
90 const RangeFieldType alpha = residual / qdoth;//<p,Bp>/<q,Aq>
91 assert( !std::isnan( alpha ) );
92 x.axpy( alpha, q );
93
94 if( preconditioner )
95 {
96 p.axpy( -alpha, h );//r_k
97 (*preconditioner)( p, s); //B*r_k
98
99 prevResidual = residual;//<rk-1,B*rk-1>
100 residual = p.scalarProductDofs( s );//<rk,B*rk>
101 }
102 else
103 {
104 r.axpy( alpha, h );
105
106 prevResidual = residual;
107 residual = r.normSquaredDofs( );
108 }
109
110 if( os )
111 {
112 (*os) << "Fem::CG it: " << iterations << " : residual " << std::sqrt(residual) << std::endl;
113 }
114 }
115
116 return (iterations < maxIterations) ? iterations : -iterations;
117 }
118
119} // naemspace LinearSolver
120} // namespace Fem
121
122} // namespace Dune
123
124#endif // #ifndef DUNE_FEM_CGINVERSEOPERATOR_HH
Type traits to determine the type of reals (when working with complex numbers)
Dune namespace.
Definition: alignedallocator.hh:13
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden & Uni Heidelberg  |  generated with Hugo v0.111.3 (Apr 19, 22:44, 2025)