dune-fem  2.4.1-rc
istlmatrixadapter.hh
Go to the documentation of this file.
1 #ifndef DUNE_FEM_ISTLMATRIXADAPTER_HH
2 #define DUNE_FEM_ISTLMATRIXADAPTER_HH
3 
4 #if HAVE_DUNE_ISTL
5 
6 #include <dune/common/timer.hh>
7 
8 #include <dune/istl/operators.hh>
9 
17 
18 namespace Dune
19 {
20  namespace Fem
21  {
22  template <class MatrixImp>
23  class LagrangeParallelMatrixAdapter
24  : public AssembledLinearOperator< MatrixImp,
25  typename MatrixImp :: RowBlockVectorType,
26  typename MatrixImp :: ColBlockVectorType>
27  {
28  typedef LagrangeParallelMatrixAdapter< MatrixImp > ThisType;
29  public:
30  typedef MatrixImp MatrixType;
31  typedef Fem::PreconditionerWrapper<MatrixType> PreconditionAdapterType;
32 
33  typedef typename MatrixType :: RowDiscreteFunctionType RowDiscreteFunctionType;
34  typedef typename MatrixType :: ColDiscreteFunctionType ColumnDiscreteFunctionType;
35 
36  typedef typename RowDiscreteFunctionType :: DiscreteFunctionSpaceType RowSpaceType;
37 
38  typedef typename ColumnDiscreteFunctionType :: DiscreteFunctionSpaceType ColSpaceType;
39  typedef Fem::ParallelScalarProduct<ColumnDiscreteFunctionType> ParallelScalarProductType;
40 
41  typedef typename RowDiscreteFunctionType :: DofStorageType X;
42  typedef typename ColumnDiscreteFunctionType :: DofStorageType Y;
43 
45  typedef MatrixType matrix_type;
46  typedef X domain_type;
47  typedef Y range_type;
48  typedef typename X::field_type field_type;
49 
51  enum { category=SolverCategory::sequential };
52 
53  public:
55  LagrangeParallelMatrixAdapter ( const LagrangeParallelMatrixAdapter &org )
56  : matrix_( org.matrix_ ),
57  rowSpace_( org.rowSpace_ ),
58  colSpace_( org.colSpace_ ),
59  scp_( colSpace_ ),
60  preconditioner_( org.preconditioner_ ),
61  averageCommTime_( org.averageCommTime_ )
62  {}
63 
66  LagrangeParallelMatrixAdapter ( MatrixType &matrix,
67  const RowSpaceType &rowSpace,
68  const ColSpaceType &colSpace,
69  const PreconditionAdapterType& precon )
70  : matrix_( matrix ),
71  rowSpace_( rowSpace ),
72  colSpace_( colSpace ),
73  scp_( colSpace_ ),
74  preconditioner_( precon ),
75  averageCommTime_( 0 )
76  {}
77 
79  double averageCommTime() const
80  {
81  return averageCommTime_ ;
82  }
83 
85  PreconditionAdapterType &preconditionAdapter()
86  {
87  return preconditioner_;
88  }
90  const PreconditionAdapterType &preconditionAdapter() const
91  {
92  return preconditioner_;
93  }
94 
96  ParallelScalarProductType &scp()
97  {
98  return scp_;
99  }
100 
102  virtual void apply ( const X &x, Y &y ) const
103  {
104  matrix_.mv( x, y );
105  communicate( y );
106  }
107 
109  virtual void applyscaleadd ( field_type alpha, const X &x, Y &y) const
110  {
111  if( rowSpace_.grid().comm().size() <= 1 )
112  {
113  matrix_.usmv(alpha,x,y);
114  communicate( y );
115  }
116  else
117  {
118  Y tmp( y.size() );
119  apply( x, tmp );
120  y.axpy( alpha, tmp );
121  }
122  }
123 
124  virtual double residuum(const Y& rhs, X& x) const
125  {
126  Y tmp( rhs );
127 
128  this->apply(x,tmp);
129  tmp -= rhs;
130 
131  // return global sum of residuum
132  return scp_.norm(tmp);
133  }
134 
136  virtual const MatrixType& getmat () const
137  {
138  return matrix_;
139  }
140 
141  protected:
142  void communicate( Y &y ) const
143  {
144  if( rowSpace_.grid().comm().size() <= 1 )
145  return;
146 
147  Dune::Timer commTime;
148  ColumnDiscreteFunctionType tmp( "LagrangeParallelMatrixAdapter::communicate", colSpace_, y );
149  colSpace_.communicate( tmp );
150  averageCommTime_ += commTime.elapsed();
151  }
152 
153  MatrixType &matrix_;
154  const RowSpaceType &rowSpace_;
155  const ColSpaceType &colSpace_;
156  mutable ParallelScalarProductType scp_;
157  PreconditionAdapterType preconditioner_;
158  mutable double averageCommTime_;
159  };
160 
161  template <class MatrixImp>
162  class DGParallelMatrixAdapter
163  : public AssembledLinearOperator< MatrixImp,
164  typename MatrixImp :: RowBlockVectorType,
165  typename MatrixImp :: ColBlockVectorType>
166  {
167  typedef DGParallelMatrixAdapter< MatrixImp > ThisType ;
168  public:
169  typedef MatrixImp MatrixType;
170  typedef Fem::PreconditionerWrapper<MatrixType> PreconditionAdapterType;
171 
172  typedef typename MatrixType :: RowDiscreteFunctionType RowDiscreteFunctionType;
173  typedef typename MatrixType :: ColDiscreteFunctionType ColumnDiscreteFunctionType;
174 
175  typedef typename RowDiscreteFunctionType :: DiscreteFunctionSpaceType RowSpaceType;
176 
177  typedef typename ColumnDiscreteFunctionType :: DiscreteFunctionSpaceType ColSpaceType;
178  typedef Fem::ParallelScalarProduct<ColumnDiscreteFunctionType> ParallelScalarProductType;
179 
180  typedef typename RowDiscreteFunctionType :: DofStorageType X;
181  typedef typename ColumnDiscreteFunctionType :: DofStorageType Y;
182 
184  typedef MatrixType matrix_type;
185  typedef X domain_type;
186  typedef Y range_type;
187  typedef typename X::field_type field_type;
188 
190  enum { category=SolverCategory::sequential };
191 
192  protected:
193  MatrixType& matrix_;
194  const RowSpaceType& rowSpace_;
195  const ColSpaceType& colSpace_;
196 
197  ParallelScalarProductType scp_;
198 
199  PreconditionAdapterType preconditioner_;
200  mutable double averageCommTime_;
201 
202  public:
204  DGParallelMatrixAdapter (const DGParallelMatrixAdapter& org)
205  : matrix_(org.matrix_)
206  , rowSpace_(org.rowSpace_)
207  , colSpace_(org.colSpace_)
208  , scp_(colSpace_)
209  , preconditioner_(org.preconditioner_)
210  , averageCommTime_( org.averageCommTime_ )
211  {}
212 
214  DGParallelMatrixAdapter (MatrixType& A,
215  const RowSpaceType& rowSpace,
216  const ColSpaceType& colSpace,
217  const PreconditionAdapterType& precon )
218  : matrix_(A)
219  , rowSpace_(rowSpace)
220  , colSpace_(colSpace)
221  , scp_(colSpace_)
222  , preconditioner_( precon )
223  , averageCommTime_( 0.0 )
224  {}
225 
227  double averageCommTime() const
228  {
229  return averageCommTime_ ;
230  }
231 
233  PreconditionAdapterType& preconditionAdapter() { return preconditioner_; }
234  const PreconditionAdapterType& preconditionAdapter() const { return preconditioner_; }
235 
237  ParallelScalarProductType& scp() { return scp_; }
238 
240  virtual void apply (const X& x, Y& y) const
241  {
242  // exchange data first
243  communicate( x );
244 
245  // apply vector to matrix
246  matrix_.mv(x,y);
247 
248  // delete non-interior
249  scp_.deleteNonInterior( y );
250  }
251 
253  virtual void applyscaleadd (field_type alpha, const X& x, Y& y) const
254  {
255  // exchange data first
256  communicate( x );
257 
258  // apply matrix
259  matrix_.usmv(alpha,x,y);
260 
261  // delete non-interior
262  scp_.deleteNonInterior( y );
263  }
264 
265  virtual double residuum(const Y& rhs, X& x) const
266  {
267  // exchange data
268  communicate( x );
269 
270  typedef typename ParallelScalarProductType :: SlaveDofsType SlaveDofsType;
271  const SlaveDofsType& slaveDofs = scp_.slaveDofs();
272 
273  typedef typename Y :: block_type LittleBlockVectorType;
274  LittleBlockVectorType tmp;
275  double res = 0.0;
276 
277  // counter for rows
278  int i = 0;
279  const int slaveSize = slaveDofs.size();
280  for(int slave = 0; slave<slaveSize; ++slave)
281  {
282  const int nextSlave = slaveDofs[slave];
283  for(; i<nextSlave; ++i)
284  {
285  tmp = 0;
286  // get row
287  typedef typename MatrixType :: row_type row_type;
288 
289  const row_type& row = matrix_[i];
290  // multiply with row
291  typedef typename MatrixType :: ConstColIterator ConstColIterator;
292  ConstColIterator endj = row.end();
293  for (ConstColIterator j = row.begin(); j!=endj; ++j)
294  {
295  (*j).umv(x[j.index()], tmp);
296  }
297 
298  // substract right hand side
299  tmp -= rhs[i];
300 
301  // add scalar product
302  res += tmp.two_norm2();
303  }
304  ++i;
305  }
306 
307  res = rowSpace_.grid().comm().sum( res );
308  // return global sum of residuum
309  return std::sqrt( res );
310  }
311 
313  virtual const MatrixType& getmat () const
314  {
315  return matrix_;
316  }
317  protected:
318  void communicate(const X& x) const
319  {
320  if( rowSpace_.grid().comm().size() <= 1 ) return ;
321 
322  Dune::Timer commTime;
323 
324  // create temporary discretet function object
325  RowDiscreteFunctionType tmp ("DGParallelMatrixAdapter::communicate",
326  rowSpace_, x );
327 
328  // exchange data by copying
329  rowSpace_.communicate( tmp );
330 
331  // accumulate communication time
332  averageCommTime_ += commTime.elapsed();
333  }
334  };
335 
336  template <class MatrixImp,class Space>
337  struct ISTLParallelMatrixAdapter;
338 
339  template< class MatrixImp,
340  class FunctionSpace, class GridPart, int polOrder, template< class > class Storage>
341  struct ISTLParallelMatrixAdapter< MatrixImp, LagrangeDiscreteFunctionSpace< FunctionSpace,GridPart,polOrder,Storage> >
342  {
343  typedef LagrangeParallelMatrixAdapter<MatrixImp> Type;
344  };
345 
346  template< class MatrixImp,
347  class FunctionSpace, class GridPart, int polOrder, template< class > class Storage>
348  struct ISTLParallelMatrixAdapter< MatrixImp, PAdaptiveLagrangeSpace< FunctionSpace,GridPart,polOrder,Storage> >
349  {
350  typedef LagrangeParallelMatrixAdapter<MatrixImp> Type;
351  };
352  template< class MatrixImp,
353  class FunctionSpace, class GridPart, int polOrder, template< class > class Storage>
354  struct ISTLParallelMatrixAdapter< MatrixImp, DiscontinuousGalerkinSpace< FunctionSpace,GridPart,polOrder,Storage> >
355  {
356  typedef DGParallelMatrixAdapter<MatrixImp> Type ;
357  };
358  template< class MatrixImp,
359  class FunctionSpace, class GridPart, int polOrder, template< class > class Storage>
360  struct ISTLParallelMatrixAdapter< MatrixImp, LagrangeDiscontinuousGalerkinSpace< FunctionSpace,GridPart,polOrder,Storage> >
361  {
362  typedef DGParallelMatrixAdapter<MatrixImp> Type ;
363  };
364  template< class MatrixImp,
365  class FunctionSpace, class GridPart, int polOrder, template< class > class Storage>
366  struct ISTLParallelMatrixAdapter< MatrixImp, LegendreDiscontinuousGalerkinSpace< FunctionSpace,GridPart,polOrder,Storage> >
367  {
368  typedef DGParallelMatrixAdapter<MatrixImp> Type ;
369  };
370  template< class MatrixImp,
371  class FunctionSpace, class GridPart, int polOrder, template< class > class Storage>
372  struct ISTLParallelMatrixAdapter< MatrixImp, HierarchicLegendreDiscontinuousGalerkinSpace< FunctionSpace,GridPart,polOrder,Storage> >
373  {
374  typedef DGParallelMatrixAdapter<MatrixImp> Type ;
375  };
376  template< class MatrixImp,
377  class FunctionSpace, class GridPart, int polOrder, template< class > class Storage>
378  struct ISTLParallelMatrixAdapter< MatrixImp, PAdaptiveDGSpace< FunctionSpace,GridPart,polOrder,Storage> >
379  {
380  typedef DGParallelMatrixAdapter<MatrixImp> Type ;
381  };
382 
383  template<class MatrixImp, class Space1, class Space2>
384  struct ISTLParallelMatrixAdapter<MatrixImp, CombinedDiscreteFunctionSpace<Space1, Space2> >
385  {
386  typedef LagrangeParallelMatrixAdapter<MatrixImp> Type;
387  };
388  template< class MatrixImp, class ... DiscreteFunctionSpaces >
389  struct ISTLParallelMatrixAdapter<MatrixImp, TupleDiscreteFunctionSpace< DiscreteFunctionSpaces ... > >
390  {
391  typedef LagrangeParallelMatrixAdapter<MatrixImp> Type;
392  };
393 
394 
395  } // end namespace Fem
396 } // end namespace Dune
397 
398 #endif // #if HAVE_DUNE_ISTL
399 
400 #endif // #ifndef DUNE_FEM_ISTLMATRIXADAPTER_HH
401 
static double sqrt(const Double &v)
Definition: double.hh:870
Definition: coordinate.hh:4