dune-istl  3.0-git
umfpack.hh
Go to the documentation of this file.
1 // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2 // vi: set et ts=4 sw=2 sts=2:
3 #ifndef DUNE_ISTL_UMFPACK_HH
4 #define DUNE_ISTL_UMFPACK_HH
5 
6 #if HAVE_SUITESPARSE_UMFPACK || defined DOXYGEN
7 
8 #include<complex>
9 #include<type_traits>
10 
11 #include<umfpack.h>
12 
13 #include<dune/common/exceptions.hh>
14 #include<dune/common/fmatrix.hh>
15 #include<dune/common/fvector.hh>
16 #include<dune/common/unused.hh>
18 #include<dune/istl/solvers.hh>
20 
21 #include"colcompmatrix.hh"
22 
23 
24 namespace Dune {
36  // FORWARD DECLARATIONS
37  template<class M, class T, class TM, class TD, class TA>
38  class SeqOverlappingSchwarz;
39 
40  template<class T, bool tag>
41  struct SeqOverlappingSchwarzAssemblerHelper;
42 
48  template<class Matrix>
49  class UMFPack
50  {};
51 
52  // wrapper class for C-Function Calls in the backend. Choose the right function namespace
53  // depending on the template parameter used.
54  template<typename T>
56  {};
57 
58  template<>
59  struct UMFPackMethodChooser<double>
60  {
61  template<typename... A>
62  static void defaults(A... args)
63  {
64  umfpack_di_defaults(args...);
65  }
66  template<typename... A>
67  static void free_numeric(A... args)
68  {
69  umfpack_di_free_numeric(args...);
70  }
71  template<typename... A>
72  static void free_symbolic(A... args)
73  {
74  umfpack_di_free_symbolic(args...);
75  }
76  template<typename... A>
77  static int load_numeric(A... args)
78  {
79  return umfpack_di_load_numeric(args...);
80  }
81  template<typename... A>
82  static void numeric(A... args)
83  {
84  umfpack_di_numeric(args...);
85  }
86  template<typename... A>
87  static void report_info(A... args)
88  {
89  umfpack_di_report_info(args...);
90  }
91  template<typename... A>
92  static void report_status(A... args)
93  {
94  umfpack_di_report_status(args...);
95  }
96  template<typename... A>
97  static int save_numeric(A... args)
98  {
99  return umfpack_di_save_numeric(args...);
100  }
101  template<typename... A>
102  static void solve(A... args)
103  {
104  umfpack_di_solve(args...);
105  }
106  template<typename... A>
107  static void symbolic(A... args)
108  {
109  umfpack_di_symbolic(args...);
110  }
111  };
112 
113  template<>
114  struct UMFPackMethodChooser<std::complex<double> >
115  {
116  template<typename... A>
117  static void defaults(A... args)
118  {
119  umfpack_zi_defaults(args...);
120  }
121  template<typename... A>
122  static void free_numeric(A... args)
123  {
124  umfpack_zi_free_numeric(args...);
125  }
126  template<typename... A>
127  static void free_symbolic(A... args)
128  {
129  umfpack_zi_free_symbolic(args...);
130  }
131  template<typename... A>
132  static int load_numeric(A... args)
133  {
134  return umfpack_zi_load_numeric(args...);
135  }
136  template<typename... A>
137  static void numeric(const int* cs, const int* ri, const double* val, A... args)
138  {
139  umfpack_zi_numeric(cs,ri,val,NULL,args...);
140  }
141  template<typename... A>
142  static void report_info(A... args)
143  {
144  umfpack_zi_report_info(args...);
145  }
146  template<typename... A>
147  static void report_status(A... args)
148  {
149  umfpack_zi_report_status(args...);
150  }
151  template<typename... A>
152  static int save_numeric(A... args)
153  {
154  return umfpack_zi_save_numeric(args...);
155  }
156  template<typename... A>
157  static void solve(int m, const int* cs, const int* ri, std::complex<double>* val, double* x, const double* b,A... args)
158  {
159  const double* cval = reinterpret_cast<const double*>(val);
160  umfpack_zi_solve(m,cs,ri,cval,NULL,x,NULL,b,NULL,args...);
161  }
162  template<typename... A>
163  static void symbolic(int m, int n, const int* cs, const int* ri, const double* val, A... args)
164  {
165  umfpack_zi_symbolic(m,n,cs,ri,val,NULL,args...);
166  }
167  };
168 
182  template<typename T, typename A, int n, int m>
183  class UMFPack<BCRSMatrix<FieldMatrix<T,n,m>,A > >
184  : public InverseOperator<
185  BlockVector<FieldVector<T,m>,
186  typename A::template rebind<FieldVector<T,m> >::other>,
187  BlockVector<FieldVector<T,n>,
188  typename A::template rebind<FieldVector<T,n> >::other> >
189  {
190  public:
193  typedef Dune::BCRSMatrix<FieldMatrix<T,n,m>,A> matrix_type;
199  typedef Dune::BlockVector<
200  FieldVector<T,m>,
201  typename A::template rebind<FieldVector<T,m> >::other> domain_type;
203  typedef Dune::BlockVector<
204  FieldVector<T,n>,
205  typename A::template rebind<FieldVector<T,n> >::other> range_type;
206 
215  UMFPack(const Matrix& matrix, int verbose=0) : matrixIsLoaded_(false)
216  {
217  //check whether T is a supported type
218  static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
219  "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
220  Caller::defaults(UMF_Control);
221  setVerbosity(verbose);
222  setMatrix(matrix);
223  }
224 
233  UMFPack(const Matrix& matrix, int verbose, bool) : matrixIsLoaded_(false)
234  {
235  //check whether T is a supported type
236  static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
237  "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
238  Caller::defaults(UMF_Control);
239  setVerbosity(verbose);
240  setMatrix(matrix);
241  }
242 
245  UMFPack() : matrixIsLoaded_(false), verbosity_(0)
246  {
247  //check whether T is a supported type
248  static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
249  "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
250  Caller::defaults(UMF_Control);
251  }
252 
263  UMFPack(const Matrix& mat_, const char* file, int verbose=0)
264  {
265  //check whether T is a supported type
266  static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
267  "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
268  Caller::defaults(UMF_Control);
269  setVerbosity(verbose);
270  int errcode = Caller::load_numeric(&UMF_Numeric, const_cast<char*>(file));
271  if ((errcode == UMFPACK_ERROR_out_of_memory) || (errcode == UMFPACK_ERROR_file_IO))
272  {
273  matrixIsLoaded_ = false;
274  setMatrix(mat_);
275  saveDecomposition(file);
276  }
277  else
278  {
279  matrixIsLoaded_ = true;
280  std::cout << "UMFPack decomposition successfully loaded from " << file << std::endl;
281  }
282  }
283 
290  UMFPack(const char* file, int verbose=0)
291  {
292  //check whether T is a supported type
293  static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
294  "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
295  Caller::defaults(UMF_Control);
296  int errcode = Caller::load_numeric(&UMF_Numeric, const_cast<char*>(file));
297  if (errcode == UMFPACK_ERROR_out_of_memory)
298  DUNE_THROW(Dune::Exception, "ran out of memory while loading UMFPack decomposition");
299  if (errcode == UMFPACK_ERROR_file_IO)
300  DUNE_THROW(Dune::Exception, "IO error while loading UMFPack decomposition");
301  matrixIsLoaded_ = true;
302  std::cout << "UMFPack decomposition successfully loaded from " << file << std::endl;
303  setVerbosity(verbose);
304  }
305 
306  virtual ~UMFPack()
307  {
308  if ((umfpackMatrix_.N() + umfpackMatrix_.M() > 0) || matrixIsLoaded_)
309  free();
310  }
311 
315  virtual void apply(domain_type& x, range_type& b, InverseOperatorResult& res)
316  {
317  if (umfpackMatrix_.N() != b.dim())
318  DUNE_THROW(Dune::ISTLError, "Size of right-hand-side vector b does not match the number of matrix rows!");
319  if (umfpackMatrix_.M() != x.dim())
320  DUNE_THROW(Dune::ISTLError, "Size of solution vector x does not match the number of matrix columns!");
321 
322  double UMF_Apply_Info[UMFPACK_INFO];
323  Caller::solve(UMFPACK_A,
324  umfpackMatrix_.getColStart(),
325  umfpackMatrix_.getRowIndex(),
326  umfpackMatrix_.getValues(),
327  reinterpret_cast<double*>(&x[0]),
328  reinterpret_cast<double*>(&b[0]),
329  UMF_Numeric,
330  UMF_Control,
331  UMF_Apply_Info);
332 
333  //this is a direct solver
334  res.iterations = 1;
335  res.converged = true;
336  res.elapsed = UMF_Apply_Info[UMFPACK_SOLVE_WALLTIME];
337 
338  printOnApply(UMF_Apply_Info);
339  }
340 
344  virtual void apply (domain_type& x, range_type& b, double reduction, InverseOperatorResult& res)
345  {
346  DUNE_UNUSED_PARAMETER(reduction);
347  apply(x,b,res);
348  }
349 
355  void apply(T* x, T* b)
356  {
357  double UMF_Apply_Info[UMFPACK_INFO];
358  Caller::solve(UMFPACK_A,
359  umfpackMatrix_.getColStart(),
360  umfpackMatrix_.getRowIndex(),
361  umfpackMatrix_.getValues(),
362  x,
363  b,
364  UMF_Numeric,
365  UMF_Control,
366  UMF_Apply_Info);
367  printOnApply(UMF_Apply_Info);
368  }
369 
381  void setOption(unsigned int option, double value)
382  {
383  if (option >= UMFPACK_CONTROL)
384  DUNE_THROW(RangeError, "Requested non-existing UMFPack option");
385 
386  UMF_Control[option] = value;
387  }
388 
392  void saveDecomposition(const char* file)
393  {
394  int errcode = Caller::save_numeric(UMF_Numeric, const_cast<char*>(file));
395  if (errcode != UMFPACK_OK)
396  DUNE_THROW(Dune::Exception,"IO ERROR while trying to save UMFPack decomposition");
397  }
398 
400  void setMatrix(const Matrix& matrix)
401  {
402  if ((umfpackMatrix_.N() + umfpackMatrix_.M() > 0) || matrixIsLoaded_)
403  free();
404  umfpackMatrix_ = matrix;
405  decompose();
406  }
407 
408  template<class S>
409  void setSubMatrix(const Matrix& _mat, const S& rowIndexSet)
410  {
411  if ((umfpackMatrix_.N() + umfpackMatrix_.M() > 0) || matrixIsLoaded_)
412  free();
413  umfpackMatrix_.setMatrix(_mat,rowIndexSet);
414  decompose();
415  }
416 
424  void setVerbosity(int v)
425  {
426  verbosity_ = v;
427  // set the verbosity level in UMFPack
428  if (verbosity_ == 0)
429  UMF_Control[UMFPACK_PRL] = 1;
430  if (verbosity_ == 1)
431  UMF_Control[UMFPACK_PRL] = 2;
432  if (verbosity_ == 2)
433  UMF_Control[UMFPACK_PRL] = 4;
434  }
435 
441  {
442  return UMF_Numeric;
443  }
444 
449  UMFPackMatrix& getInternalMatrix()
450  {
451  return umfpackMatrix_;
452  }
453 
458  void free()
459  {
460  if (!matrixIsLoaded_)
461  {
462  Caller::free_symbolic(&UMF_Symbolic);
463  umfpackMatrix_.free();
464  }
465  Caller::free_numeric(&UMF_Numeric);
466  matrixIsLoaded_ = false;
467  }
468 
469  const char* name() { return "UMFPACK"; }
470 
471  private:
472  typedef typename Dune::UMFPackMethodChooser<T> Caller;
473 
474  template<class M,class X, class TM, class TD, class T1>
475  friend class SeqOverlappingSchwarz;
476  friend struct SeqOverlappingSchwarzAssemblerHelper<UMFPack<Matrix>,true>;
477 
479  void decompose()
480  {
481  double UMF_Decomposition_Info[UMFPACK_INFO];
482  Caller::symbolic(static_cast<int>(umfpackMatrix_.N()),
483  static_cast<int>(umfpackMatrix_.N()),
484  umfpackMatrix_.getColStart(),
485  umfpackMatrix_.getRowIndex(),
486  reinterpret_cast<double*>(umfpackMatrix_.getValues()),
487  &UMF_Symbolic,
488  UMF_Control,
489  UMF_Decomposition_Info);
490  Caller::numeric(umfpackMatrix_.getColStart(),
491  umfpackMatrix_.getRowIndex(),
492  reinterpret_cast<double*>(umfpackMatrix_.getValues()),
493  UMF_Symbolic,
494  &UMF_Numeric,
495  UMF_Control,
496  UMF_Decomposition_Info);
497  Caller::report_status(UMF_Control,UMF_Decomposition_Info[UMFPACK_STATUS]);
498  if (verbosity_ == 1)
499  {
500  std::cout << "[UMFPack Decomposition]" << std::endl;
501  std::cout << "Wallclock Time taken: " << UMF_Decomposition_Info[UMFPACK_NUMERIC_WALLTIME] << " (CPU Time: " << UMF_Decomposition_Info[UMFPACK_NUMERIC_TIME] << ")" << std::endl;
502  std::cout << "Flops taken: " << UMF_Decomposition_Info[UMFPACK_FLOPS] << std::endl;
503  std::cout << "Peak Memory Usage: " << UMF_Decomposition_Info[UMFPACK_PEAK_MEMORY]*UMF_Decomposition_Info[UMFPACK_SIZE_OF_UNIT] << " bytes" << std::endl;
504  std::cout << "Condition number estimate: " << 1./UMF_Decomposition_Info[UMFPACK_RCOND] << std::endl;
505  std::cout << "Numbers of non-zeroes in decomposition: L: " << UMF_Decomposition_Info[UMFPACK_LNZ] << " U: " << UMF_Decomposition_Info[UMFPACK_UNZ] << std::endl;
506  }
507  if (verbosity_ == 2)
508  {
509  Caller::report_info(UMF_Control,UMF_Decomposition_Info);
510  }
511  }
512 
513  void printOnApply(double* UMF_Info)
514  {
515  Caller::report_status(UMF_Control,UMF_Info[UMFPACK_STATUS]);
516  if (verbosity_ > 0)
517  {
518  std::cout << "[UMFPack Solve]" << std::endl;
519  std::cout << "Wallclock Time: " << UMF_Info[UMFPACK_SOLVE_WALLTIME] << " (CPU Time: " << UMF_Info[UMFPACK_SOLVE_TIME] << ")" << std::endl;
520  std::cout << "Flops Taken: " << UMF_Info[UMFPACK_SOLVE_FLOPS] << std::endl;
521  std::cout << "Iterative Refinement steps taken: " << UMF_Info[UMFPACK_IR_TAKEN] << std::endl;
522  std::cout << "Error Estimate: " << UMF_Info[UMFPACK_OMEGA1] << " resp. " << UMF_Info[UMFPACK_OMEGA2] << std::endl;
523  }
524  }
525 
526  UMFPackMatrix umfpackMatrix_;
527  bool matrixIsLoaded_;
528  int verbosity_;
529  void *UMF_Symbolic;
530  void *UMF_Numeric;
531  double UMF_Control[UMFPACK_CONTROL];
532  };
533 
534  template<typename T, typename A, int n, int m>
536  {
537  enum { value=true};
538  };
539 
540  template<typename T, typename A, int n, int m>
542  {
543  enum { value = true };
544  };
545 }
546 
547 #endif // HAVE_SUITESPARSE_UMFPACK
548 
549 #endif //DUNE_ISTL_UMFPACK_HH
Definition: solvertype.hh:27
UMFPack()
default constructor
Definition: umfpack.hh:245
virtual void apply(domain_type &x, range_type &b, InverseOperatorResult &res)
Apply inverse operator,.
Definition: umfpack.hh:315
Use the UMFPack package to directly solve linear systems – empty default class.
Definition: umfpack.hh:49
UMFPack(const Matrix &matrix, int verbose=0)
Construct a solver object from a BCRSMatrix.
Definition: umfpack.hh:215
void apply(T *x, T *b)
additional apply method with c-arrays in analogy to superlu
Definition: umfpack.hh:355
UMFPackMatrix & getInternalMatrix()
Return the column compress matrix from UMFPack.
Definition: umfpack.hh:449
static int load_numeric(A... args)
Definition: umfpack.hh:77
Statistics about the application of an inverse operator.
Definition: solver.hh:31
Dune::BlockVector< FieldVector< T, m >, typename A::template rebind< FieldVector< T, m > >::other > domain_type
The type of the domain of the solver.
Definition: umfpack.hh:201
static void free_numeric(A... args)
Definition: umfpack.hh:67
static void free_symbolic(A... args)
Definition: umfpack.hh:72
double elapsed
Elapsed time in seconds.
Definition: solver.hh:62
Implementations of the inverse operator interface.
Inititializer for the ColCompMatrix as needed by OverlappingSchwarz.
Definition: colcompmatrix.hh:153
static int save_numeric(A... args)
Definition: umfpack.hh:152
static void report_info(A... args)
Definition: umfpack.hh:142
static void report_status(A... args)
Definition: umfpack.hh:92
size_type dim() const
dimension of the vector space
Definition: bvector.hh:279
Templates characterizing the type of a solver.
static void numeric(const int *cs, const int *ri, const double *val, A... args)
Definition: umfpack.hh:137
virtual void apply(domain_type &x, range_type &b, double reduction, InverseOperatorResult &res)
apply inverse operator, with given convergence criteria.
Definition: umfpack.hh:344
void saveDecomposition(const char *file)
saves a decomposition to a file
Definition: umfpack.hh:392
Definition: colcompmatrix.hh:160
Definition: solvertype.hh:13
const char * name()
Definition: umfpack.hh:469
Dune::BlockVector< FieldVector< T, n >, typename A::template rebind< FieldVector< T, n > >::other > range_type
The type of the range of the solver.
Definition: umfpack.hh:205
static void free_numeric(A... args)
Definition: umfpack.hh:122
virtual ~UMFPack()
Definition: umfpack.hh:306
Definition: matrixutils.hh:25
UMFPack(const Matrix &mat_, const char *file, int verbose=0)
Try loading a decomposition from file and do a decomposition if unsuccessful.
Definition: umfpack.hh:263
A sparse block matrix with compressed row storage.
Definition: bcrsmatrix.hh:411
static void numeric(A... args)
Definition: umfpack.hh:82
int iterations
Number of iterations.
Definition: solver.hh:50
Dune::ColCompMatrix< Matrix > UMFPackMatrix
The corresponding SuperLU Matrix type.
Definition: umfpack.hh:195
Definition: basearray.hh:19
static void report_info(A... args)
Definition: umfpack.hh:87
static void defaults(A... args)
Definition: umfpack.hh:62
bool converged
True if convergence criterion has been met.
Definition: solver.hh:56
void setMatrix(const Matrix &matrix)
Initialize data from given matrix.
Definition: umfpack.hh:400
static void solve(A... args)
Definition: umfpack.hh:102
Implementation of the BCRSMatrix class.
A vector of blocks with memory management.
Definition: bvector.hh:308
void * getFactorization()
Return the matrix factorization.
Definition: umfpack.hh:440
static int save_numeric(A... args)
Definition: umfpack.hh:97
Definition: umfpack.hh:55
void setSubMatrix(const Matrix &_mat, const S &rowIndexSet)
Definition: umfpack.hh:409
static void solve(int m, const int *cs, const int *ri, std::complex< double > *val, double *x, const double *b, A... args)
Definition: umfpack.hh:157
void free()
free allocated space.
Definition: umfpack.hh:458
static void report_status(A... args)
Definition: umfpack.hh:147
static void symbolic(int m, int n, const int *cs, const int *ri, const double *val, A... args)
Definition: umfpack.hh:163
Abstract base class for all solvers.
Definition: solver.hh:79
UMFPack(const char *file, int verbose=0)
try loading a decomposition from file
Definition: umfpack.hh:290
void setVerbosity(int v)
sets the verbosity level for the UMFPack solver
Definition: umfpack.hh:424
UMFPack(const Matrix &matrix, int verbose, bool)
Constructor for compatibility with SuperLU standard constructor.
Definition: umfpack.hh:233
derive error class from the base class in common
Definition: istlexception.hh:16
Sequential overlapping Schwarz preconditioner.
Definition: colcompmatrix.hh:157
static int load_numeric(A... args)
Definition: umfpack.hh:132
STL namespace.
static void symbolic(A... args)
Definition: umfpack.hh:107
static void defaults(A... args)
Definition: umfpack.hh:117
static void free_symbolic(A... args)
Definition: umfpack.hh:127
void setOption(unsigned int option, double value)
Set UMFPack-specific options.
Definition: umfpack.hh:381