dune-fem  2.4.1-rc
oemsolver/oemsolver.hh
Go to the documentation of this file.
1 #ifndef DUNE_FEM_OEMSOLVER_HH
2 #define DUNE_FEM_OEMSOLVER_HH
3 
4 //- system includes
5 #include <limits>
6 #include <utility>
7 
8 //- Dune includes
11 
12 //- local includes
13 #include "preconditioning.hh"
14 
15 #include <dune/fem/solver/pardg.hh>
16 
17 // include BLAS implementation
18 #include "cblas.h"
19 
20 
21 namespace Dune
22 {
23 
24  namespace Fem
25  {
26 
31  struct OEMMatrix
32  {
33  virtual ~OEMMatrix () {}
34 
40  virtual void multOEM ( const double *u, double *w ) const = 0;
41 
47  virtual double ddotOEM ( const double *u, const double *v ) const = 0;
48 
49  // SparseRowMatrixObject does not implement this method
50 #if 0
51 
56  virtual void precondition ( const double *u, double *w ) const = 0;
57 #endif
58  };
59 
60  } // end namespace Fem
61 
62  } // end namespace Dune
63 
64 
65 
66  namespace OEMSolver
67  {
68 
70  //
71  // Operator Interface to use linear solvers from pardg
72  //
74  template <class OperatorImp>
77  : public PARDG::Function
78 #endif
79  {
80  const OperatorImp & op_;
81  int size_;
82  public:
83  SolverInterfaceImpl(const OperatorImp & op, int size = 0)
84  : op_(op), size_(size)
85  {}
86 
87  void setSize( int size ) { size_ = size; }
88 
89  void operator () (const double *arg, double * dest, int i = 0 )
90  {
91  op_.multOEM(arg,dest);
92  }
93 
94  void mult(const double *arg, double * dest) const
95  {
96  op_.multOEM(arg,dest);
97  }
98 
99  int dim_of_argument(int i = 0) const
100  {
101  assert( i == 0 );
102  return size_;
103  }
104  int dim_of_value(int i = 0) const
105  {
106  assert( i == 0 );
107  return size_;
108  }
109  };
110 
112  //
113  // Preconditioner Interface to use linear solvers from pardg
114  //
116  template <class PreconditionerImp>
119  : public PARDG::Function
120 #endif
121  {
122  const PreconditionerImp& pre_;
123  int size_;
124  public:
125  PreconditionerImpl(const PreconditionerImp& pre, int size = 0)
126  : pre_(pre), size_(size)
127  {}
128 
129  void setSize( int size ) { size_ = size; }
131  void operator () (const double *arg, double * dest, int i = 0 )
132  {
133  pre_.precondition(arg,dest);
134  }
135 
136  void mult(const double *arg, double * dest) const
137  {
138  pre_.precondition(arg,dest);
139  }
141  int dim_of_argument(int i = 0) const
142  {
143  assert( i == 0 );
144  return size_;
145  }
146  int dim_of_value(int i = 0) const
147  {
148  assert( i == 0 );
149  return size_;
150  }
151  };
152 
153  // use cblas implementations
154  using namespace DuneCBlas;
155 
156  using DuneCBlas :: daxpy;
157  using DuneCBlas :: dcopy;
158  using DuneCBlas :: ddot;
159  using DuneCBlas :: dnrm2;
160  using DuneCBlas :: dscal;
161 
164  template <class MatrixImp, class VectorType>
165  void mult(const MatrixImp & m, const VectorType * x, VectorType * ret)
166  {
167  // call multOEM of the matrix
168  m.multOEM(x,ret);
169  }
170 
172  template <class Matrix , class PC_Matrix , bool >
173  struct Mult
174  {
175  static inline double ddot( const Matrix& A,
176  const double *x,
177  const double *y)
178  {
179  return A.ddotOEM(x,y);
180  }
181 
182  typedef bool mult_t(const Matrix &A,
183  const PC_Matrix & C,
184  const double *arg,
185  double *dest ,
186  double * tmp);
187 
188  static void back_solve(const int size,
189  const PC_Matrix & C, double* solution, double* tmp)
190  {
191  assert( tmp );
192  if( C.rightPrecondition() )
193  {
194  C.precondition(solution,tmp);
195  // copy modified solution
196  std::memcpy(solution,tmp, size * sizeof(double));
197  }
198  }
200  static bool mult_pc (const Matrix &A, const PC_Matrix & C,
201  const double *arg, double *dest , double * tmp)
202  {
203  assert( tmp );
204 
205  bool rightPreCon = C.rightPrecondition();
206  // check type of preconditioning
207  if( rightPreCon )
208  {
209  // call precondition of Matrix PC
210  C.precondition(arg,tmp);
211 
212  // call mult of Matrix A
213  mult(A,tmp,dest);
214  }
215  else
216  {
217  // call mult of Matrix A
218  mult(A,arg,tmp);
219 
220  // call precondition of Matrix PC
221  C.precondition(tmp,dest);
222  }
223  return rightPreCon ;
224  }
225  };
226 
228  template <class Matrix>
229  struct Mult<Matrix,Matrix,false>
230  {
231  static inline double ddot( const Matrix& A,
232  const double *x,
233  const double *y)
234  {
235  return A.ddotOEM(x,y);
236  }
237 
238  typedef bool mult_t(const Matrix &A,
239  const Matrix &C,
240  const double *arg,
241  double *dest ,
242  double * tmp);
243 
244  static void back_solve(const int size,
245  const Matrix & C, double* solution, double* tmp)
246  {
247  // do nothing here
248  }
249 
250  static bool mult_pc(const Matrix &A, const Matrix & C, const double *arg ,
251  double *dest , double * tmp)
252  {
253  // tmp has to be 0
254  assert( tmp == 0 );
255  // C is just a fake
256  assert( &A == &C );
257 
258  // call mult of Matrix A
259  mult(A,arg,dest);
260 
261  return true;
262  }
263  };
264 
265 #define USE_MEMPROVIDER
266 #include "bicgstab.h"
267 #include "cghs.h"
268 #include "gmres.h"
269 #include "bicgsq.h"
270 #undef USE_MEMPROVIDER
271 
272 
276  {
277  // size of vectors
278  const int size_;
279 
280  // indices of external values
281  std::vector<int> indices_;
282  public:
283  // use with care, not sure that working correctly already
284  template <class SolverOperatorImp>
285  FakeConditioner(int size, SolverOperatorImp& op) : size_(size)
286  {
287  assert( size_ > 0 );
288 
289  double * diag = new double [size_];
290  double * tmp = new double [size_];
291 
292  assert( diag );
293  assert( tmp );
294  for(int i=0; i<size_; ++i) tmp[i] = i;
295  op(tmp,diag);
296 
297  int newSize = (int) 0.25 * size_;
298  indices_.reserve( newSize );
299  indices_.resize( 0 );
300  // now diag contains only non-zeros for all internal entries
301  // these are set to 1.0 to be the id mapping
302  for(int i=0; i<size_; ++i)
303  {
304  if( ! (std::abs (diag[i]) > 0.0) )
305  {
306  indices_.push_back( i );
307  }
308  }
309 
310  delete [] diag;
311  delete [] tmp;
312  }
313 
314  bool rightPrecondition() const { return false; }
315 
317  void precondition(const double * arg, double * dest) const
318  {
319  multOEM(arg,dest);
320  }
321 
323  void multOEM(const double * arg, double * dest) const
324  {
325  std::memcpy( dest, arg , size_ * sizeof(double) );
326 
327  const int s = indices_.size();
328  for(int i=0; i<s; ++i)
329  {
330  dest[indices_[i]] = 0.0;
331  }
332  }
333  };
334 
335  } // end namespace OEMSolver
336 
337  namespace Dune
338  {
339  namespace Fem
340  {
352  template <class DiscreteFunctionType, class OpType >
353  class OEMCGOp : public Operator<
354  DiscreteFunctionType,DiscreteFunctionType>
355  {
356  public:
357  typedef OpType OperatorType;
358 
359  private:
360  // no const reference, we make const later
361  OperatorType &op_;
362  typename DiscreteFunctionType::RangeFieldType epsilon_;
363  int maxIter_;
364  bool verbose_ ;
365  mutable int iterations_;
366 
367  typedef std::pair < int , double > ReturnValueType;
368 
369  template <class OperatorImp, bool hasPreconditioning>
370  struct SolverCaller
371  {
372  template <class DiscreteFunctionImp>
373  static ReturnValueType call(OperatorImp & op,
374  const DiscreteFunctionImp & arg,
375  DiscreteFunctionImp & dest,
376  double eps, int maxIter, bool verbose)
377  {
378  // use communication class of grid
379  // see dune-common/common/collectivecommunication.hh
380  // for interface
381  int size = arg.space().size();
382 
383  if(op.hasPreconditionMatrix())
384  {
385  return OEMSolver::cghs(arg.space().gridPart().comm(),
386  size,op.systemMatrix(),op.preconditionMatrix(),
387  arg.leakPointer(),dest.leakPointer(),eps,maxIter,verbose);
388  }
389  else
390  {
391  return OEMSolver::cghs(arg.space().gridPart().comm(),
392  size,op.systemMatrix(),
393  arg.leakPointer(),dest.leakPointer(),eps,maxIter,verbose);
394  }
395  }
396  };
397 
399  template <class OperatorImp>
400  struct SolverCaller<OperatorImp,false>
401  {
402  template <class DiscreteFunctionImp>
403  static ReturnValueType call(OperatorImp & op,
404  const DiscreteFunctionImp & arg,
405  DiscreteFunctionImp & dest,
406  double eps, int maxIter, bool verbose)
407  {
408  // use communication class of grid
409  // see dune-common/common/collectivecommunication.hh
410  // for interface
411  int size = arg.space().size();
412  return OEMSolver::cghs(arg.space().gridPart().comm(),
413  size,op.systemMatrix(),
414  arg.leakPointer(),dest.leakPointer(),eps,maxIter,verbose);
415  }
416  };
417 
418  public:
419 
427  OEMCGOp ( OperatorType &op,
428  double redEps, double absLimit, int maxIter, bool verbose,
429  const ParameterReader &parameter = Parameter::container() )
430  : op_(op),
431  epsilon_( absLimit ),
432  maxIter_( maxIter ),
433  verbose_( verbose ),
434  iterations_( 0 )
435  {}
436 
437  OEMCGOp ( OperatorType &op,
438  double redEps, double absLimit, int maxIter,
439  const ParameterReader &parameter = Parameter::container() )
440  : op_( op ),
441  epsilon_( absLimit ),
442  maxIter_( maxIter ),
443  verbose_( parameter.getValue< bool >( "fem.solver.verbose", false ) ),
444  iterations_( 0 )
445  {}
446 
447  OEMCGOp ( OperatorType &op,
448  double redEps, double absLimit,
449  const ParameterReader &parameter = Parameter::container() )
450  : op_( op ),
451  epsilon_( absLimit ),
452  maxIter_( std::numeric_limits< int >::max() ),
453  verbose_( parameter.getValue< bool >( "fem.solver.verbose", false ) ),
454  iterations_( 0 )
455  {}
456 
457  void prepare (const DiscreteFunctionType& Arg, DiscreteFunctionType& Dest) const
458  {
459  }
460 
461  void finalize () const
462  {
463  }
464 
465  int iterations () const
466  {
467  return iterations_;
468  }
469 
474  void apply( const DiscreteFunctionType& arg, DiscreteFunctionType& dest ) const
475  {
476  // prepare operator
477  prepare ( arg, dest );
478 
479  ReturnValueType val =
480  SolverCaller<OperatorType,
481  // check wheter operator has precondition methods
482  // to enable preconditioning derive your operator from
483  // OEMSolver::PreconditionInterface
484  Conversion<OperatorType, OEMSolver::PreconditionInterface > ::exists >::
485  // call solver, see above
486  call(op_,arg,dest,epsilon_,maxIter_,verbose_);
487 
488  iterations_ = val.first;
489 
490  if(arg.space().gridPart().comm().rank() == 0)
491  {
492  std::cout << "OEM-CG: " << val.first << " iterations! Error: " << val.second << "\n";
493  }
494 
495  // finalize operator
496  finalize ();
497  }
498 
503  void operator ()( const DiscreteFunctionType& arg, DiscreteFunctionType& dest ) const
504  {
505  apply(arg,dest);
506  }
507  };
508 
510  template <class DiscreteFunctionType, class OpType>
511  class OEMBICGSTABOp : public Operator<
512  DiscreteFunctionType,DiscreteFunctionType>
513  {
514  public:
515  typedef OpType OperatorType;
516 
517  private:
518  // no const reference, we make const later
519  OperatorType &op_;
520  typename DiscreteFunctionType::RangeFieldType epsilon_;
521  int maxIter_;
522  bool verbose_ ;
523  mutable int iterations_;
524 
525  typedef std::pair < int , double > ReturnValueType;
526 
527  template <class OperatorImp, bool hasPreconditioning>
528  struct SolverCaller
529  {
530  template <class DiscreteFunctionImp>
531  static ReturnValueType call(OperatorImp & op,
532  const DiscreteFunctionImp & arg,
533  DiscreteFunctionImp & dest,
534  double eps, int maxIter, bool verbose)
535  {
536  int size = arg.space().size();
537  if(op.hasPreconditionMatrix())
538  {
539  return OEMSolver::bicgstab(arg.space().gridPart().comm(),
540  size,op.systemMatrix(),op.preconditionMatrix(),
541  arg.leakPointer(),dest.leakPointer(),eps,maxIter,verbose);
542  }
543  else
544  {
545  return OEMSolver::bicgstab(arg.space().gridPart().comm(),
546  size,op.systemMatrix(),
547  arg.leakPointer(),dest.leakPointer(),eps,maxIter,verbose);
548  }
549  }
550  };
551 
553  template <class OperatorImp>
554  struct SolverCaller<OperatorImp,false>
555  {
556  template <class DiscreteFunctionImp>
557  static ReturnValueType call(OperatorImp & op,
558  const DiscreteFunctionImp & arg,
559  DiscreteFunctionImp & dest,
560  double eps, int maxIter, bool verbose)
561  {
562  int size = arg.space().size();
563  return OEMSolver::bicgstab(arg.space().gridPart().comm(),
564  size,op.systemMatrix(),
565  arg.leakPointer(),dest.leakPointer(),eps,maxIter,verbose);
566  }
567  };
568 
569  public:
577  OEMBICGSTABOp ( OperatorType& op,
578  double redEps, double absLimit, int maxIter, bool verbose,
579  const ParameterReader &parameter = Parameter::container() )
580  : op_(op),
581  epsilon_( absLimit ),
582  maxIter_( maxIter ),
583  verbose_( verbose ),
584  iterations_( 0 )
585  {
586  }
587 
588  OEMBICGSTABOp ( OperatorType &op,
589  double redEps, double absLimit, int maxIter,
590  const ParameterReader &parameter = Parameter::container() )
591  : op_( op ),
592  epsilon_( absLimit ),
593  maxIter_( maxIter ),
594  verbose_( parameter.getValue< bool >( "fem.solver.verbose", false ) ),
595  iterations_( 0 )
596  {
597  }
598 
599  OEMBICGSTABOp ( OperatorType &op,
600  double redEps, double absLimit,
601  const ParameterReader &parameter = Parameter::container() )
602  : op_( op ),
603  epsilon_( absLimit ),
604  maxIter_( std::numeric_limits< int >::max() ),
605  verbose_( parameter.getValue< bool >( "fem.solver.verbose", false ) ),
606  iterations_( 0 )
607  {
608  }
609  void prepare (const DiscreteFunctionType& Arg, DiscreteFunctionType& Dest) const
610  {
611  }
612 
613  void finalize () const
614  {
615  }
616 
617  int iterations () const
618  {
619  return iterations_;
620  }
621 
622 
627  void apply( const DiscreteFunctionType& arg, DiscreteFunctionType& dest ) const
628  {
629  // prepare operator
630  prepare ( arg, dest );
631 
632  ReturnValueType val =
633  SolverCaller<OperatorType,
634  // check wheter operator has precondition methods
635  // to enable preconditioning derive your operator from
636  // OEMSolver::PreconditionInterface
637  Conversion<OperatorType, OEMSolver::PreconditionInterface > ::exists >::
638  // call solver, see above
639  call(op_,arg,dest,epsilon_,maxIter_,verbose_);
640 
641  iterations_ = val.first;
642 
643  if(arg.space().gridPart().comm().rank() == 0)
644  {
645  std::cout << "OEM-BICGstab: " << val.first << " iterations! Error: " << val.second << "\n";
646  }
647 
648  // finalize operator
649  finalize ();
650  }
651 
656  void operator ()( const DiscreteFunctionType& arg, DiscreteFunctionType& dest ) const
657  {
658  apply(arg,dest);
659  }
660 
661  };
662 
664  // BICG SQ scheme
666 
667  template <class DiscreteFunctionType, class OpType>
668  class OEMBICGSQOp : public Operator<
669  DiscreteFunctionType,DiscreteFunctionType>
670  {
671  public:
672  typedef OpType OperatorType;
673 
674  private:
675  // no const reference, we make const later
676  OperatorType &op_;
677  typename DiscreteFunctionType::RangeFieldType epsilon_;
678  int maxIter_;
679  bool verbose_ ;
680  mutable int iterations_;
681 
682  public:
690  OEMBICGSQOp ( OperatorType &op,
691  double redEps, double absLimit, int maxIter, bool verbose,
692  const ParameterReader &parameter = Parameter::container() )
693  : op_(op),
694  epsilon_( absLimit ),
695  maxIter_( maxIter ),
696  verbose_( verbose ),
697  iterations_( 0 )
698  {
699  }
700 
701  OEMBICGSQOp ( OperatorType &op,
702  double redEps, double absLimit,
703  const ParameterReader &parameter = Parameter::container() )
704  : op_( op ),
705  epsilon_( absLimit ),
706  maxIter_( std::numeric_limits< int >::max() ),
707  verbose_( parameter.getValue< bool >( "fem.solver.verbose", false ) ),
708  iterations_( 0 )
709  {
710  }
711 
712  OEMBICGSQOp ( OperatorType &op,
713  double redEps, double absLimit, int maxIter,
714  const ParameterReader &parameter = Parameter::container() )
715  : op_( op ),
716  epsilon_( absLimit ),
717  maxIter_( maxIter ),
718  verbose_( parameter.getValue< bool >( "fem.solver.verbose", false ) ),
719  iterations_( 0 )
720  {
721  }
722 
723  void prepare (const DiscreteFunctionType& Arg, DiscreteFunctionType& Dest) const
724  {
725  }
726 
727  void finalize () const
728  {
729  }
730 
731  int iterations () const
732  {
733  return iterations_;
734  }
735 
740  void apply( const DiscreteFunctionType& arg, DiscreteFunctionType& dest ) const
741  {
742  // prepare operator
743  prepare ( arg, dest );
744 
745  int size = arg.space().size();
746 
747  int iter = OEMSolver::bicgsq(size,op_.systemMatrix(),
748  arg.leakPointer(),dest.leakPointer(),epsilon_,maxIter_,verbose_);
749 
750  iterations_ = iter;
751 
752  std::cout << "OEM-BICGGsq: " << iter << " iterations!\n";
753  // finalize operator
754  finalize ();
755  }
756 
761  void operator ()( const DiscreteFunctionType& arg, DiscreteFunctionType& dest ) const
762  {
763  apply(arg,dest);
764  }
765 
766  };
767 
768 
770  template< class DiscreteFunctionType, class Op >
772  : public Operator< DiscreteFunctionType, DiscreteFunctionType >
773  {
775 
776  public:
777  typedef Op OperatorType;
778 
779  private:
780  // type of internal projector if no preconditioner given
782 
783  typedef std::pair < int , double > ReturnValueType;
784 
785  template <class OperatorImp, bool hasPreconditioning>
786  struct SolverCaller
787  {
788  template <class DiscreteFunctionImp>
789  static ReturnValueType call(OperatorImp & op,
790  const DiscreteFunctionImp & arg,
791  DiscreteFunctionImp & dest,
792  int inner, double eps, int maxIter, bool verbose)
793  {
794  int size = arg.space().size();
795  if(op.hasPreconditionMatrix())
796  {
797  return OEMSolver::gmres(arg.space().gridPart().comm(),
798  inner,size,op.systemMatrix(),op.preconditionMatrix(),
799  arg.leakPointer(),dest.leakPointer(),eps,maxIter,verbose);
800  }
801  // in parallel case we need special treatment, if no preconditoner exist
802  else if( arg.space().gridPart().comm().size() > 1 )
803  {
804  typedef typename OperatorImp::SystemMatrixType SystemMatrixType;
805  OEMSolver::SolverInterfaceImpl< SystemMatrixType > opSolve( op.systemMatrix() );
806  FakeConditionerType preConditioner( size, opSolve );
807  return OEMSolver::gmres(arg.space().gridPart().comm(),
808  inner,size,op.systemMatrix(),preConditioner,
809  arg.leakPointer(),dest.leakPointer(),eps,maxIter,verbose);
810  }
811  else
812  {
813  return OEMSolver::gmres(arg.space().gridPart().comm(),
814  inner,size,op.systemMatrix(),
815  arg.leakPointer(),dest.leakPointer(),eps,maxIter,verbose);
816  }
817  }
818  };
819 
820  // without any preconditioning
821  template <class OperatorImp>
822  struct SolverCaller<OperatorImp,false>
823  {
824  template <class DiscreteFunctionImp>
825  static ReturnValueType call(OperatorImp & op,
826  const DiscreteFunctionImp & arg,
827  DiscreteFunctionImp & dest,
828  int inner, double eps, int maxIter, bool verbose)
829  {
830  int size = arg.space().size();
831  if( arg.space().gridPart().comm().size() > 1 )
832  {
833  typedef typename OperatorImp::MatrixType MatrixType;
834  OEMSolver::SolverInterfaceImpl< MatrixType > opSolve( op.matrix() );
835  FakeConditionerType preConditioner( size, opSolve );
836  return OEMSolver::gmres(arg.space().gridPart().comm(),
837  inner,size,op.matrix(),preConditioner,
838  arg.leakPointer(),dest.leakPointer(),eps,maxIter,verbose);
839  }
840  else
841  {
842  return OEMSolver::gmres(arg.space().gridPart().comm(),
843  inner,size,op.matrix(),
844  arg.leakPointer(),dest.leakPointer(),eps,maxIter,verbose);
845  }
846  }
847  };
848 
849  public:
857  OEMGMRESOp ( OperatorType &op,
858  double redEps, double absLimit, int maxIter, bool verbose,
859  const ParameterReader &parameter = Parameter::container() )
860  : op_( op ),
861  epsilon_( absLimit ),
862  maxIter_( maxIter ),
863  restart_( parameter.getValue< int >( "oemsolver.gmres.restart", 20 ) ),
864  verbose_( verbose ),
865  iterations_( 0 )
866  {}
867 
868  OEMGMRESOp ( OperatorType &op,
869  double redEps, double absLimit, int maxIter,
870  const ParameterReader &parameter = Parameter::container() )
871  : op_( op ),
872  epsilon_( absLimit ),
873  maxIter_( maxIter ),
874  restart_( parameter.getValue< int >( "oemsolver.gmres.restart", 20 ) ),
875  verbose_( parameter.getValue< bool >( "fem.solver.verbose", false ) ),
876  iterations_( 0 )
877  {}
878 
879  OEMGMRESOp ( OperatorType &op,
880  double redEps, double absLimit,
881  const ParameterReader &parameter = Parameter::container() )
882  : op_( op ),
883  epsilon_( absLimit ),
884  maxIter_( std::numeric_limits< int >::max() ),
885  restart_( parameter.getValue< int >( "oemsolver.gmres.restart", 20 ) ),
886  verbose_( parameter.getValue< bool >( "fem.solver.verbose", false ) ),
887  iterations_( 0 )
888  {}
889 
890  void prepare (const DiscreteFunctionType& Arg, DiscreteFunctionType& Dest) const
891  {
892  }
893 
894  void finalize () const
895  {
896  }
897 
898  int iterations() const
899  {
900  return iterations_;
901  }
902 
907  void apply( const DiscreteFunctionType& arg, DiscreteFunctionType& dest ) const
908  {
909  // prepare operator
910  prepare ( arg, dest );
911 
912  int size = arg.space().size();
913  int inner = std::min( size, restart_ );
914 
915  ReturnValueType val =
916  SolverCaller<OperatorType,
917  // check wheter operator has precondition methods
918  // to enable preconditioning derive your operator from
919  // OEMSolver::PreconditionInterface
920  Conversion<OperatorType, OEMSolver::PreconditionInterface > ::exists >::
921  // call solver, see above
922  call(op_,arg,dest,inner,epsilon_,maxIter_,verbose_);
923 
924  iterations_ = val.first;
925 
926  if( arg.space().gridPart().comm().rank() == 0 && verbose_ )
927  {
928  std::cout << "OEM-GMRES: " << val.first << " iterations! Error: " << val.second << "\n";
929  }
930 
931  // finalize operator
932  finalize ();
933  }
934 
939  void operator ()( const DiscreteFunctionType& arg, DiscreteFunctionType& dest ) const
940  {
941  apply(arg,dest);
942  }
943 
944  private:
945  // no const reference, we make const later
946  OperatorType &op_;
947  typename DiscreteFunctionType::RangeFieldType epsilon_;
948  int maxIter_;
949  int restart_;
950  bool verbose_ ;
951  mutable int iterations_;
952  };
953 
957 #ifdef USE_PARDG_ODE_SOLVER
958  //
960  // GMRES Version of Dennis code
961  //
963  // \brief GMRES implementation from Dennis D.
964  template <class DiscreteFunctionType, class OperatorImp>
965  class GMRESOp : public Operator<
966  DiscreteFunctionType,DiscreteFunctionType>
967  {
968  public:
969  typedef OperatorImp OperatorType;
970  private:
971  typedef OEMSolver :: FakeConditioner FakeConditioner;
972 
973  template <class SolverType, bool hasPreconditioning>
974  struct SolverCaller
975  {
976  template< class OperatorImpA, class PreConMatrix, class DiscreteFunctionImp >
977  static void solve(SolverType & solver,
978  OperatorImpA & op,
979  const PreConMatrix & pm,
980  const DiscreteFunctionImp & arg,
981  DiscreteFunctionImp & dest)
982  {
983  int size = arg.space().size();
984 
986 
987  // in parallel runs we need fake pre conditioner to
988  // project vectors onto interior
989  if(op.hasPreconditionMatrix())
990  {
992  solver.set_preconditioner(pre);
993 
994  // note argument and destination are toggled
995  solver.solve(opSolve, dest.leakPointer() , arg.leakPointer() );
996 
997  solver.unset_preconditioner();
998  }
999  else
1000  {
1001  // note argument and destination are toggled
1002  solver.solve(opSolve, dest.leakPointer() , arg.leakPointer() );
1003  }
1004  }
1005 
1006  template <class OperatorImpA, class DiscreteFunctionImp>
1007  static void call(SolverType & solver,
1008  OperatorImpA & op,
1009  const DiscreteFunctionImp & arg,
1010  DiscreteFunctionImp & dest)
1011  {
1012  solve(solver,op,op.preconditionMatrix(),arg,dest);
1013  }
1014  };
1015 
1016  // without any preconditioning
1017  template <class SolverType>
1018  struct SolverCaller<SolverType,false>
1019  {
1020  template <class OperatorImpA, class DiscreteFunctionImp>
1021  static void call(SolverType & solver,
1022  OperatorImpA & op,
1023  const DiscreteFunctionImp & arg,
1024  DiscreteFunctionImp & dest)
1025  {
1026  int size = arg.space().size();
1028 
1029  // in parallel runs we need fake pre conditioner to
1030  // project vectors onto interior
1031  if(arg.space().gridPart().comm().size() > 1)
1032  {
1033  FakeConditioner fake(size,opSolve);
1035  solver.set_preconditioner(pre);
1036 
1037  // note argument and destination are toggled
1038  solver.solve(opSolve, dest.leakPointer() , arg.leakPointer() );
1039  solver.unset_preconditioner();
1040  }
1041  else
1042  {
1043  // note argument and destination are toggled
1044  solver.solve(opSolve, dest.leakPointer() , arg.leakPointer() );
1045  }
1046  }
1047  };
1048 
1049  // solver
1050  typedef PARDG::GMRES SolverType;
1051  mutable SolverType solver_;
1052 
1053  // wrapper to fit interface of FGMRES operator
1054  OperatorType &op_;
1055 
1056  typename DiscreteFunctionType::RangeFieldType epsilon_;
1057  int maxIter_;
1058  bool verbose_ ;
1059  mutable int iterations_;
1060 
1061  public:
1062  GMRESOp( OperatorType & op , double redEps , double absLimit , int maxIter , bool verbose,
1063  const ParameterReader &parameter = Parameter::container() )
1064  : solver_(PARDG::Communicator::instance(), parameter.getValue< int >( "oemsolver.gmres.restart", 20 ) )
1065  , op_(op) , epsilon_ ( absLimit )
1066  , maxIter_ (maxIter ) , verbose_ ( verbose )
1067  , iterations_ ( 0 )
1068  {
1069  }
1070 
1071  void prepare (const DiscreteFunctionType& Arg, DiscreteFunctionType& Dest) const
1072  {
1073  }
1074 
1075  void finalize () const
1076  {
1077  }
1078 
1079  int iterations () const
1080  {
1081  return iterations_;
1082  }
1083 
1088  void apply( const DiscreteFunctionType& arg, DiscreteFunctionType& dest ) const
1089  {
1090  // prepare operator
1091  prepare ( arg, dest );
1092 
1093  solver_.set_tolerance(epsilon_);
1094  solver_.set_max_number_of_iterations( maxIter_ );
1095 
1096  if(verbose_)
1097  {
1098  solver_.IterativeSolver::set_output(std::cout);
1099  solver_.DynamicalObject::set_output(std::cout);
1100  }
1101 
1102  SolverCaller<SolverType,
1103  // check wheter operator has precondition methods
1104  // to enable preconditioning derive your operator from
1105  // OEMSolver::PreconditionInterface
1106  Conversion<OperatorType, OEMSolver::PreconditionInterface > ::exists >::
1107  // call solver, see above
1108  call(solver_,op_.systemMatrix(),arg,dest);
1109 
1110  iterations_ = solver_.number_of_iterations();
1111 
1112  // finalize operator
1113  finalize ();
1114  }
1115 
1120  void operator ()( const DiscreteFunctionType& arg, DiscreteFunctionType& dest ) const
1121  {
1122  apply(arg,dest);
1123  }
1124  };
1125 
1126  template <class DiscreteFunctionType, class OperatorImp>
1127  class FGMRESOp : public Operator<
1128  DiscreteFunctionType,DiscreteFunctionType>
1129  {
1130  public:
1131  typedef OperatorImp OperatorType;
1132  private:
1133  typedef OEMSolver :: FakeConditioner FakeConditionerType;
1134 
1135  template <class SolverType, bool hasPreconditioning>
1136  struct SolverCaller
1137  {
1138  template <class OperatorImpA, class PreConMatrix, class DiscreteFunctionImp>
1139  static void solve(SolverType & solver,
1140  OperatorImpA & op,
1141  const PreConMatrix & pm,
1142  const DiscreteFunctionImp & arg,
1143  DiscreteFunctionImp & dest)
1144  {
1145  int size = arg.space().size();
1148  solver.set_preconditioner(pre);
1149 
1150  // note argument and destination are toggled
1151  solver.solve(opSolve, dest.leakPointer() , arg.leakPointer() );
1152  solver.unset_preconditioner();
1153  }
1154 
1155  template <class OperatorImpA, class DiscreteFunctionImp>
1156  static void call(SolverType & solver,
1157  OperatorImpA & op,
1158  const DiscreteFunctionImp & arg,
1159  DiscreteFunctionImp & dest)
1160  {
1161  if(op.hasPreconditionMatrix() )
1162  {
1163  solve(solver,op.systemMatrix(),op.preconditionMatrix(),arg,dest);
1164  }
1165  else
1166  {
1167  SolverCaller<SolverType,false>::call(solver,op,arg,dest);
1168  }
1169  }
1170  };
1171 
1172  // without any preconditioning
1173  template <class SolverType>
1174  struct SolverCaller<SolverType,false>
1175  {
1176  template <class OperatorImpA, class DiscreteFunctionImp>
1177  static void solve(SolverType & solver,
1178  OperatorImpA & op,
1179  const DiscreteFunctionImp & arg,
1180  DiscreteFunctionImp & dest)
1181  {
1182  int size = arg.space().size();
1184  FakeConditionerType fake(size,opSolve);
1185  SolverCaller<SolverType,true>::solve(solver,op,fake,arg,dest);
1186  }
1187 
1188  template <class OperatorImpA, class DiscreteFunctionImp>
1189  static void call(SolverType & solver,
1190  OperatorImpA & op,
1191  const DiscreteFunctionImp & arg,
1192  DiscreteFunctionImp & dest)
1193  {
1194  // not working yet
1195  assert( false );
1196  solve(solver,op.systemMatrix(),arg,dest);
1197  }
1198  };
1199 
1200  // solver
1201  typedef PARDG::FGMRES SolverType;
1202  mutable SolverType solver_;
1203 
1204  // wrapper to fit interface of FGMRES operator
1205  OperatorType &op_;
1206 
1207  typename DiscreteFunctionType::RangeFieldType epsilon_;
1208  int maxIter_;
1209  bool verbose_ ;
1210  mutable int iterations_;
1211 
1212  public:
1213  FGMRESOp( OperatorType & op , double redEps , double absLimit , int maxIter , bool verbose,
1214  const ParameterReader &parameter = Parameter::container() )
1215  : solver_(PARDG::Communicator::instance(), parameter.getValue< int >( "oemsolver.gmres.restart", 20 ) )
1216  , op_(op) , epsilon_ ( absLimit )
1217  , maxIter_ (maxIter ) , verbose_ ( verbose )
1218  , iterations_( 0 )
1219  {
1220  }
1221 
1222  void prepare (const DiscreteFunctionType& Arg, DiscreteFunctionType& Dest) const
1223  {
1224  }
1225 
1226  void finalize () const
1227  {
1228  }
1229 
1230  int iterations () const
1231  {
1232  return iterations_;
1233  }
1234 
1236  void apply( const DiscreteFunctionType& arg, DiscreteFunctionType& dest ) const
1237  {
1238  // prepare operator
1239  prepare ( arg, dest );
1240 
1241  solver_.set_tolerance(epsilon_);
1242  solver_.set_max_number_of_iterations( maxIter_ );
1243 
1244  if(verbose_)
1245  {
1246  solver_.IterativeSolver::set_output(std::cout);
1247  solver_.DynamicalObject::set_output(std::cout);
1248  }
1249 
1250  SolverCaller<SolverType,
1251  // check wheter operator has precondition methods
1252  // to enable preconditioning derive your operator from
1253  // OEMSolver::PreconditionInterface
1254  Conversion<OperatorType, OEMSolver::PreconditionInterface > ::exists >::
1255  // call solver, see above
1256  call(solver_,op_,arg,dest);
1257 
1258  iterations_ = solver_.number_of_iterations();
1259 
1260  // finalize operator
1261  finalize ();
1262  }
1263 
1265  void operator ()( const DiscreteFunctionType& arg, DiscreteFunctionType& dest ) const
1266  {
1267  apply(arg,dest);
1268  }
1269 
1270  };
1271 
1273  //
1274  // BICGstab Version of Dennis code
1275  //
1277  /*
1278  \interface
1279  \brief BICG-stab implementation from Dennis D.
1280  */
1281  template <class DiscreteFunctionType, class OperatorImp>
1282  class BICGSTABOp : public Operator<
1283  DiscreteFunctionType,DiscreteFunctionType>
1284  {
1285  public:
1286  typedef OperatorImp OperatorType;
1287  private:
1288  template <class SolverType, bool hasPreconditioning>
1289  struct SolverCaller
1290  {
1291  template <class OperatorImpA, class PreConMatrix, class DiscreteFunctionImp>
1292  static void solve(SolverType & solver,
1293  OperatorImpA & op,
1294  const PreConMatrix & pm,
1295  const DiscreteFunctionImp & arg,
1296  DiscreteFunctionImp & dest)
1297  {
1298  int size = arg.space().size();
1300 
1302  solver.set_preconditioner(pre);
1303 
1304  // note argument and destination are toggled
1305  solver.solve(opSolve, dest.leakPointer() , arg.leakPointer() );
1306  solver.unset_preconditioner();
1307  }
1308 
1309  template <class OperatorImpA, class DiscreteFunctionImp>
1310  static void call(SolverType & solver,
1311  OperatorImpA & op,
1312  const DiscreteFunctionImp & arg,
1313  DiscreteFunctionImp & dest)
1314  {
1315  if(op.hasPreconditionMatrix())
1316  {
1317  solve(solver,op.systemMatrix(),op.preconditionMatrix(),arg,dest);
1318  }
1319  else
1320  {
1321  SolverCaller<SolverType,false>::call(solver,op,arg,dest);
1322  }
1323  }
1324  };
1325 
1326  // without any preconditioning
1327  template <class SolverType>
1328  struct SolverCaller<SolverType,false>
1329  {
1330  template <class OperatorImpA, class DiscreteFunctionImp>
1331  static void solve(SolverType & solver,
1332  OperatorImpA & op,
1333  const DiscreteFunctionImp & arg,
1334  DiscreteFunctionImp & dest)
1335  {
1336  int size = arg.space().size();
1338 
1339  // note argument and destination are toggled
1340  solver.solve(opSolve, dest.leakPointer() , arg.leakPointer() );
1341  }
1342  template <class OperatorImpA, class DiscreteFunctionImp>
1343  static void call(SolverType & solver,
1344  OperatorImpA & op,
1345  const DiscreteFunctionImp & arg,
1346  DiscreteFunctionImp & dest)
1347  {
1348  solve(solver,op.systemMatrix(),arg,dest);
1349  }
1350  };
1351 
1352  // solver
1353  typedef PARDG::BICGSTAB SolverType;
1354  mutable SolverType solver_;
1355  // wrapper to fit interface of GMRES operator
1356  OperatorType &op_;
1357 
1358  typename DiscreteFunctionType::RangeFieldType epsilon_;
1359  int maxIter_;
1360  bool verbose_ ;
1361  mutable int iterations_;
1362 
1363  public:
1364  BICGSTABOp( OperatorType & op , double redEps , double absLimit , int maxIter , bool verbose ,
1365  const ParameterReader &parameter = Parameter::container() )
1366  : solver_(PARDG::Communicator::instance())
1367  , op_(op), epsilon_ ( absLimit )
1368  , maxIter_ (maxIter ) , verbose_ ( verbose )
1369  , iterations_( 0 )
1370  {
1371  }
1372 
1373  void prepare (const DiscreteFunctionType& Arg, DiscreteFunctionType& Dest) const
1374  {
1375  }
1376 
1377  void finalize () const
1378  {
1379  }
1380 
1381  int iterations () const
1382  {
1383  return iterations_;
1384  }
1385 
1387  void apply( const DiscreteFunctionType& arg, DiscreteFunctionType& dest ) const
1388  {
1389  // prepare operator
1390  prepare ( arg, dest );
1391 
1392  solver_.set_tolerance( epsilon_ );
1393  solver_.set_max_number_of_iterations( maxIter_ );
1394 
1395  if( verbose_ )
1396  {
1397  solver_.IterativeSolver::set_output(std::cout);
1398  solver_.DynamicalObject::set_output(std::cout);
1399  }
1400 
1401  SolverCaller<SolverType,
1402  // check wheter operator has precondition methods
1403  // to enable preconditioning derive your operator from
1404  // OEMSolver::PreconditionInterface
1405  Conversion<OperatorType, OEMSolver::PreconditionInterface > ::exists >::
1406  // call solver, see above
1407  call(solver_,op_,arg,dest);
1408 
1409  iterations_ = solver_.number_of_iterations();
1410 
1411  // finalize operator
1412  finalize ();
1413  }
1414 
1416  void operator ()( const DiscreteFunctionType& arg, DiscreteFunctionType& dest ) const
1417  {
1418  apply(arg,dest);
1419  }
1420 
1421  };
1422 #endif
1423 
1424  } // namespace Fem
1425 
1426 } // namespace Dune
1427 
1428 #endif //#indef DUNE_FEM_OEMSOLVER_HH
void mult(const double *arg, double *dest) const
Definition: oemsolver/oemsolver.hh:94
int iterations() const
Definition: oemsolver/oemsolver.hh:617
Op OperatorType
Definition: oemsolver/oemsolver.hh:777
static bool mult_pc(const Matrix &A, const PC_Matrix &C, const double *arg, double *dest, double *tmp)
Definition: oemsolver/oemsolver.hh:200
static constexpr T min(T a)
Definition: utility.hh:81
void apply(const DiscreteFunctionType &arg, DiscreteFunctionType &dest) const
solve the system
Definition: oemsolver/oemsolver.hh:627
void setSize(int size)
Definition: oemsolver/oemsolver.hh:87
void apply(const DiscreteFunctionType &arg, DiscreteFunctionType &dest) const
solve the system
Definition: oemsolver/oemsolver.hh:907
std::pair< int, double > bicgstab(const CommunicatorType &comm, unsigned int N, const MATRIX &A, const double *b, double *x, double eps, int maxIter, bool verbose)
Definition: oemsolver/oemsolver.hh:199
mult method when given pre conditioning matrix
Definition: oemsolver/oemsolver.hh:173
Definition: oemsolver/oemsolver.hh:75
static void back_solve(const int size, const PC_Matrix &C, double *solution, double *tmp)
Definition: oemsolver/oemsolver.hh:188
OEMBICGSQOp(OperatorType &op, double redEps, double absLimit, const ParameterReader &parameter=Parameter::container())
Definition: oemsolver/oemsolver.hh:701
bool mult_t(const Matrix &A, const Matrix &C, const double *arg, double *dest, double *tmp)
Definition: oemsolver/oemsolver.hh:238
static void back_solve(const int size, const Matrix &C, double *solution, double *tmp)
Definition: oemsolver/oemsolver.hh:244
OEM-CG scheme after Hestenes and Stiefel.
Definition: oemsolver/oemsolver.hh:353
OEMBICGSQOp(OperatorType &op, double redEps, double absLimit, int maxIter, bool verbose, const ParameterReader &parameter=Parameter::container())
constructor of OEM-BiCG-SQ
Definition: oemsolver/oemsolver.hh:690
int dim_of_argument(int i=0) const
Definition: oemsolver/oemsolver.hh:99
static double max(const Double &v, const double p)
Definition: double.hh:387
PreconditionerImpl(const PreconditionerImp &pre, int size=0)
Definition: oemsolver/oemsolver.hh:125
void setSize(int size)
Definition: oemsolver/oemsolver.hh:129
OEMCGOp(OperatorType &op, double redEps, double absLimit, int maxIter, bool verbose, const ParameterReader &parameter=Parameter::container())
constructor of OEM-CG
Definition: oemsolver/oemsolver.hh:427
OEMCGOp(OperatorType &op, double redEps, double absLimit, int maxIter, const ParameterReader &parameter=Parameter::container())
Definition: oemsolver/oemsolver.hh:437
static double ddot(const Matrix &A, const double *x, const double *y)
Definition: oemsolver/oemsolver.hh:175
OEMGMRESOp(OperatorType &op, double redEps, double absLimit, const ParameterReader &parameter=Parameter::container())
Definition: oemsolver/oemsolver.hh:879
FakeConditioner(int size, SolverOperatorImp &op)
Definition: oemsolver/oemsolver.hh:285
int dim_of_value(int i=0) const
Definition: oemsolver/oemsolver.hh:104
void multOEM(const double *arg, double *dest) const
only keep internal parts of arg
Definition: oemsolver/oemsolver.hh:323
void prepare(const DiscreteFunctionType &Arg, DiscreteFunctionType &Dest) const
Definition: oemsolver/oemsolver.hh:890
Double abs(const Double &a)
Definition: double.hh:860
abstract operator
Definition: operator.hh:25
#define USE_PARDG_ODE_SOLVER
Definition: pardg.hh:15
int dim_of_value(int i=0) const
Definition: oemsolver/oemsolver.hh:146
virtual ~OEMMatrix()
Definition: oemsolver/oemsolver.hh:33
BiCG-stab solver.
Definition: oemsolver/oemsolver.hh:511
Definition: coordinate.hh:4
std::pair< int, double > cghs(const CommunicatorType &comm, unsigned int N, const MATRIX &A, const PC_MATRIX &C, const double *b, double *x, double eps, int maxIter, bool detailed)
Definition: oemsolver/oemsolver.hh:139
OpType OperatorType
Definition: oemsolver/oemsolver.hh:515
void finalize() const
Definition: oemsolver/oemsolver.hh:613
Definition: oemsolver/oemsolver.hh:275
static bool mult_pc(const Matrix &A, const Matrix &C, const double *arg, double *dest, double *tmp)
Definition: oemsolver/oemsolver.hh:250
void mult(const double *arg, double *dest) const
Definition: oemsolver/oemsolver.hh:136
static ParameterContainer & container()
Definition: io/parameter.hh:190
STL namespace.
std::pair< int, double > gmres(const CommunicatorType &comm, int m, int n, const Matrix &A, const double *b, double *x, double eps, int maxIter, bool verbose)
Definition: oemsolver/oemsolver.hh:203
static double ddot(const Matrix &A, const double *x, const double *y)
Definition: oemsolver/oemsolver.hh:231
virtual double ddotOEM(const double *u, const double *v) const =0
evaluate scalar product
void prepare(const DiscreteFunctionType &Arg, DiscreteFunctionType &Dest) const
Definition: oemsolver/oemsolver.hh:457
bool mult_t(const Matrix &A, const PC_Matrix &C, const double *arg, double *dest, double *tmp)
Definition: oemsolver/oemsolver.hh:182
OEMBICGSTABOp(OperatorType &op, double redEps, double absLimit, const ParameterReader &parameter=Parameter::container())
Definition: oemsolver/oemsolver.hh:599
OEMBICGSQOp(OperatorType &op, double redEps, double absLimit, int maxIter, const ParameterReader &parameter=Parameter::container())
Definition: oemsolver/oemsolver.hh:712
Definition: oemsolver/oemsolver.hh:66
int iterations() const
Definition: oemsolver/oemsolver.hh:731
bool rightPrecondition() const
Definition: oemsolver/oemsolver.hh:314
interface for matrices to be used with OEM sovlers
Definition: oemsolver/oemsolver.hh:31
OEMBICGSTABOp(OperatorType &op, double redEps, double absLimit, int maxIter, bool verbose, const ParameterReader &parameter=Parameter::container())
constructor of OEM-BiCG-stab
Definition: oemsolver/oemsolver.hh:577
Definition: oemsolver/oemsolver.hh:117
int iterations() const
Definition: oemsolver/oemsolver.hh:465
void apply(const DiscreteFunctionType &arg, DiscreteFunctionType &dest) const
solve the system
Definition: oemsolver/oemsolver.hh:740
OEMCGOp(OperatorType &op, double redEps, double absLimit, const ParameterReader &parameter=Parameter::container())
Definition: oemsolver/oemsolver.hh:447
int dim_of_argument(int i=0) const
Definition: oemsolver/oemsolver.hh:141
int bicgsq(unsigned N, const MATRIX &A, const double *b, double *x, double eps, int maxIter)
Definition: oemsolver/oemsolver.hh:202
OpType OperatorType
Definition: oemsolver/oemsolver.hh:672
void finalize() const
Definition: oemsolver/oemsolver.hh:461
void prepare(const DiscreteFunctionType &Arg, DiscreteFunctionType &Dest) const
Definition: oemsolver/oemsolver.hh:723
virtual void multOEM(const double *u, double *w) const =0
evaluate matrix vector multiplication
void prepare(const DiscreteFunctionType &Arg, DiscreteFunctionType &Dest) const
Definition: oemsolver/oemsolver.hh:609
void mult(const MatrixImp &m, const VectorType *x, VectorType *ret)
Definition: oemsolver/oemsolver.hh:165
BiCG-SQ method.
Definition: oemsolver/oemsolver.hh:668
SolverInterfaceImpl(const OperatorImp &op, int size=0)
Definition: oemsolver/oemsolver.hh:83
OEMGMRESOp(OperatorType &op, double redEps, double absLimit, int maxIter, const ParameterReader &parameter=Parameter::container())
Definition: oemsolver/oemsolver.hh:868
void finalize() const
Definition: oemsolver/oemsolver.hh:894
OEMBICGSTABOp(OperatorType &op, double redEps, double absLimit, int maxIter, const ParameterReader &parameter=Parameter::container())
Definition: oemsolver/oemsolver.hh:588
OpType OperatorType
Definition: oemsolver/oemsolver.hh:357
OEMGMRESOp(OperatorType &op, double redEps, double absLimit, int maxIter, bool verbose, const ParameterReader &parameter=Parameter::container())
constructor of OEM-GMRES
Definition: oemsolver/oemsolver.hh:857
void precondition(const double *arg, double *dest) const
only keep internal parts of arg
Definition: oemsolver/oemsolver.hh:317
void apply(const DiscreteFunctionType &arg, DiscreteFunctionType &dest) const
solve the system
Definition: oemsolver/oemsolver.hh:474
GMRES solver.
Definition: oemsolver/oemsolver.hh:771
void finalize() const
Definition: oemsolver/oemsolver.hh:727
int iterations() const
Definition: oemsolver/oemsolver.hh:898