dune-fem  2.4.1-rc
preconditionerwrapper.hh
Go to the documentation of this file.
1 #ifndef DUNE_FEM_PRECONDITIONERWRAPPER_HH
2 #define DUNE_FEM_PRECONDITIONERWRAPPER_HH
3 
4 #include <memory>
5 #include <dune/common/shared_ptr.hh>
6 
7 // standard diagonal preconditioner
9 
10 
11 #if HAVE_DUNE_ISTL
12 #include <dune/istl/operators.hh>
13 #include <dune/istl/preconditioners.hh>
14 
15 #include <dune/istl/paamg/amg.hh>
16 #include <dune/istl/paamg/pinfo.hh>
17 #endif
18 
19 namespace Dune
20 {
21 
22  namespace Fem
23  {
24 
25 #if HAVE_DUNE_ISTL
26  template<class M, class X, class Y, int l=1>
27  class FemSeqILU0 : public SeqILU0<M,X,Y,l>
28  {
29  typedef SeqILU0<M,X,Y,l> BaseType ;
30  public:
31  typedef typename X::field_type field_type;
32  FemSeqILU0 (const M& A, int iter, field_type w)
33  : BaseType( A, w )
34  {}
35  FemSeqILU0 (const M& A, field_type w)
36  : BaseType( A, w )
37  {}
38  };
39 
40 
41  template< class MatrixObject, class X, class Y >
42  class FemDiagonalPreconditioner : public Preconditioner<X,Y>
43  {
44  public:
45  typedef typename MatrixObject :: ColumnDiscreteFunctionType DiscreteFunctionType ;
46  // select preconditioner for assembled operators
47  typedef DiagonalPreconditionerBase< DiscreteFunctionType, MatrixObject, true > PreconditionerType ;
48 
49  protected:
50  PreconditionerType diagonalPrecon_;
51 
52  public:
53  typedef typename X::field_type field_type;
54  FemDiagonalPreconditioner( const MatrixObject& mObj )
55  : diagonalPrecon_( mObj )
56  {}
57 
59  virtual void pre (X& x, Y& b) {}
60 
62  virtual void apply (X& v, const Y& d)
63  {
64  diagonalPrecon_.applyToISTLBlockVector( d, v );
65  }
66 
68  virtual void post (X& x) {}
69  };
70 
71 
75  template<class MatrixImp>
76  class PreconditionerWrapper
77  : public Preconditioner<typename MatrixImp :: RowBlockVectorType,
78  typename MatrixImp :: ColBlockVectorType>
79  {
80  typedef MatrixImp MatrixType;
81  typedef typename MatrixType :: RowBlockVectorType X;
82  typedef typename MatrixType :: ColBlockVectorType Y;
83 
84  // use BCRSMatrix type because of specializations in dune-istl
85  typedef typename MatrixType :: BaseType ISTLMatrixType ;
86 
87  typedef typename MatrixType :: CollectiveCommunictionType
88  CollectiveCommunictionType;
89 
90  // matrix adapter for AMG
91  //#if HAVE_MPI
92  // typedef Dune::OverlappingSchwarzOperator<
93  // ISTLMatrixType, X, Y, CollectiveCommunictionType> OperatorType ;
94  //#else
95  typedef MatrixAdapter< ISTLMatrixType, X, Y > OperatorType;
96  //#endif
97  mutable std::shared_ptr< OperatorType > op_;
98 
99  // auto pointer to preconditioning object
100  typedef Preconditioner<X,Y> PreconditionerInterfaceType;
101  mutable std::shared_ptr<PreconditionerInterfaceType> preconder_;
102 
103  // flag whether we have preconditioning, and if yes if it is AMG
104  const int preEx_;
105 
106  template <class XImp, class YImp>
107  struct Apply
108  {
109  inline static void apply(PreconditionerInterfaceType& preconder,
110  XImp& v, const YImp& d)
111  {
112  }
113 
114  inline static void copy(XImp& v, const YImp& d)
115  {
116  }
117  };
118 
119  template <class XImp>
120  struct Apply<XImp,XImp>
121  {
122 
123  inline static void apply(PreconditionerInterfaceType& preconder,
124  XImp& v, const XImp& d)
125  {
126  preconder.apply( v ,d );
127  }
128 
129  inline static void copy(X& v, const Y& d)
130  {
131  v = d;
132  }
133  };
134  public:
136  typedef X domain_type;
138  typedef Y range_type;
140  typedef typename X::field_type field_type;
141 
142  enum {
144  category=SolverCategory::sequential };
145 
147  PreconditionerWrapper (const PreconditionerWrapper& org)
148  : op_( org.op_ )
149  , preconder_(org.preconder_)
150  , preEx_(org.preEx_)
151  {
152  }
153 
155  PreconditionerWrapper()
156  : op_()
157  , preconder_()
158  , preEx_( 0 )
159  {}
160 
162  template <class PreconditionerType>
163  PreconditionerWrapper(MatrixType & matrix,
164  int iter,
165  field_type relax,
166  const PreconditionerType* p )
167  : op_()
168  , preconder_( new PreconditionerType( matrix, iter, relax ) )
169  , preEx_( 1 )
170  {
171  }
172 
173 
176  template <class PreconditionerType>
177  PreconditionerWrapper(MatrixType & matrix,
178  PreconditionerType* p )
179  : op_()
180  , preconder_( p )
181  , preEx_( 1 )
182  {
183  }
184 
186  template <class PreconditionerType>
187  PreconditionerWrapper(MatrixType & matrix,
188  int iter,
189  field_type relax,
190  const PreconditionerType* p ,
191  const CollectiveCommunictionType& comm )
192  //#if HAVE_MPI
193  // : op_( new OperatorType( matrix, comm ) )
194  //#else
195  : op_( new OperatorType( matrix ) )
196  //#endif
197  , preconder_( createAMGPreconditioner(comm, iter, relax, p) )
198  , preEx_( 2 )
199  {
200  }
201 
203  virtual void pre (X& x, Y& b)
204  {
205  // all the sequentiel implemented Preconditioners do nothing in pre and post
206  // apply preconditioner
207  if( preEx_ > 1 )
208  {
209  //X tmp (x);
210  preconder_->pre(x,b);
211  //assert( std::abs( x.two_norm() - tmp.two_norm() ) < 1e-15);
212  }
213  }
214 
216  virtual void apply (X& v, const Y& d)
217  {
218  if( preEx_ )
219  {
220  // apply preconditioner
221  Apply<X,Y>::apply( *preconder_ , v, d );
222  }
223  else
224  {
225  Apply<X,Y>::copy( v, d );
226  }
227  }
228 
230  virtual void post (X& x)
231  {
232  // all the implemented Preconditioners do nothing in pre and post
233  // apply preconditioner
234  if( preEx_ > 1 )
235  {
236  preconder_->post(x);
237  }
238  }
239 
240  protected:
241  template <class Smoother>
242  PreconditionerInterfaceType*
243  createAMGPreconditioner(const CollectiveCommunictionType& comm,
244  int iter, field_type relax, const Smoother* )
245  {
246  typedef typename Dune::FieldTraits< field_type>::real_type real_type;
247  typedef typename std::conditional< std::is_convertible<field_type,real_type>::value,
248  Dune::Amg::FirstDiagonal, Dune::Amg::RowSum >::type Norm;
249  typedef Dune::Amg::CoarsenCriterion<
250  Dune::Amg::UnSymmetricCriterion<ISTLMatrixType, Norm > > Criterion;
251 
252  typedef typename Dune::Amg::SmootherTraits<Smoother>::Arguments SmootherArgs;
253 
254  SmootherArgs smootherArgs;
255 
256  smootherArgs.iterations = iter;
257  smootherArgs.relaxationFactor = relax ;
258 
259  int coarsenTarget=1200;
260  Criterion criterion(15,coarsenTarget);
261  criterion.setDefaultValuesIsotropic(2);
262  criterion.setAlpha(.67);
263  criterion.setBeta(1.0e-8);
264  criterion.setMaxLevel(10);
265 
266  /*
267  if( comm.size() > 1 )
268  {
269  typedef Dune::OwnerOverlapCopyCommunication<int> ParallelInformation;
270  ParallelInformation pinfo(MPI_COMM_WORLD);
271  typedef Dune::Amg::AMG<OperatorType, X, Smoother, ParallelInformation> AMG;
272  return new AMG(*op_, criterion, smootherArgs, pinfo);
273  }
274  else
275  */
276  {
277  // X == Y is needed for AMG
278  typedef Dune::Amg::AMG<OperatorType, X, Smoother> AMG;
279  return new AMG(*op_, criterion, smootherArgs);
280  //return new AMG(*op_, criterion, smootherArgs, 1, 1, 1, false);
281  }
282  }
283  };
284 #endif // end HAVE_DUNE_ISTL
285 
286  } // namespace Fem
287 
288 } // namespace Dune
289 
290 #endif // #ifndef DUNE_FEM_PRECONDITIONERWRAPPER_HH
Definition: coordinate.hh:4