Dune Core Modules (2.5.2)

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 
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>
55  struct UMFPackMethodChooser
56  {
57  static constexpr bool valid = false ;
58  };
59 
60  template<>
61  struct UMFPackMethodChooser<double>
62  {
63  static constexpr bool valid = true ;
64 
65  template<typename... A>
66  static void defaults(A... args)
67  {
68  umfpack_di_defaults(args...);
69  }
70  template<typename... A>
71  static void free_numeric(A... args)
72  {
73  umfpack_di_free_numeric(args...);
74  }
75  template<typename... A>
76  static void free_symbolic(A... args)
77  {
78  umfpack_di_free_symbolic(args...);
79  }
80  template<typename... A>
81  static int load_numeric(A... args)
82  {
83  return umfpack_di_load_numeric(args...);
84  }
85  template<typename... A>
86  static void numeric(A... args)
87  {
88  umfpack_di_numeric(args...);
89  }
90  template<typename... A>
91  static void report_info(A... args)
92  {
93  umfpack_di_report_info(args...);
94  }
95  template<typename... A>
96  static void report_status(A... args)
97  {
98  umfpack_di_report_status(args...);
99  }
100  template<typename... A>
101  static int save_numeric(A... args)
102  {
103  return umfpack_di_save_numeric(args...);
104  }
105  template<typename... A>
106  static void solve(A... args)
107  {
108  umfpack_di_solve(args...);
109  }
110  template<typename... A>
111  static void symbolic(A... args)
112  {
113  umfpack_di_symbolic(args...);
114  }
115  };
116 
117  template<>
118  struct UMFPackMethodChooser<std::complex<double> >
119  {
120  static constexpr bool valid = true ;
121 
122  template<typename... A>
123  static void defaults(A... args)
124  {
125  umfpack_zi_defaults(args...);
126  }
127  template<typename... A>
128  static void free_numeric(A... args)
129  {
130  umfpack_zi_free_numeric(args...);
131  }
132  template<typename... A>
133  static void free_symbolic(A... args)
134  {
135  umfpack_zi_free_symbolic(args...);
136  }
137  template<typename... A>
138  static int load_numeric(A... args)
139  {
140  return umfpack_zi_load_numeric(args...);
141  }
142  template<typename... A>
143  static void numeric(const int* cs, const int* ri, const double* val, A... args)
144  {
145  umfpack_zi_numeric(cs,ri,val,NULL,args...);
146  }
147  template<typename... A>
148  static void report_info(A... args)
149  {
150  umfpack_zi_report_info(args...);
151  }
152  template<typename... A>
153  static void report_status(A... args)
154  {
155  umfpack_zi_report_status(args...);
156  }
157  template<typename... A>
158  static int save_numeric(A... args)
159  {
160  return umfpack_zi_save_numeric(args...);
161  }
162  template<typename... A>
163  static void solve(int m, const int* cs, const int* ri, std::complex<double>* val, double* x, const double* b,A... args)
164  {
165  const double* cval = reinterpret_cast<const double*>(val);
166  umfpack_zi_solve(m,cs,ri,cval,NULL,x,NULL,b,NULL,args...);
167  }
168  template<typename... A>
169  static void symbolic(int m, int n, const int* cs, const int* ri, const double* val, A... args)
170  {
171  umfpack_zi_symbolic(m,n,cs,ri,val,NULL,args...);
172  }
173  };
174 
188  template<typename T, typename A, int n, int m>
189  class UMFPack<BCRSMatrix<FieldMatrix<T,n,m>,A > >
190  : public InverseOperator<
191  BlockVector<FieldVector<T,m>,
192  typename A::template rebind<FieldVector<T,m> >::other>,
193  BlockVector<FieldVector<T,n>,
194  typename A::template rebind<FieldVector<T,n> >::other> >
195  {
196  public:
205  typedef Dune::BlockVector<
207  typename A::template rebind<FieldVector<T,m> >::other> domain_type;
209  typedef Dune::BlockVector<
211  typename A::template rebind<FieldVector<T,n> >::other> range_type;
212 
221  UMFPack(const Matrix& matrix, int verbose=0) : matrixIsLoaded_(false)
222  {
223  //check whether T is a supported type
224  static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
225  "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
226  Caller::defaults(UMF_Control);
227  setVerbosity(verbose);
228  setMatrix(matrix);
229  }
230 
239  UMFPack(const Matrix& matrix, int verbose, bool) : matrixIsLoaded_(false)
240  {
241  //check whether T is a supported type
242  static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
243  "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
244  Caller::defaults(UMF_Control);
245  setVerbosity(verbose);
246  setMatrix(matrix);
247  }
248 
251  UMFPack() : matrixIsLoaded_(false), verbosity_(0)
252  {
253  //check whether T is a supported type
254  static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
255  "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
256  Caller::defaults(UMF_Control);
257  }
258 
269  UMFPack(const Matrix& mat_, const char* file, int verbose=0)
270  {
271  //check whether T is a supported type
272  static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
273  "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
274  Caller::defaults(UMF_Control);
275  setVerbosity(verbose);
276  int errcode = Caller::load_numeric(&UMF_Numeric, const_cast<char*>(file));
277  if ((errcode == UMFPACK_ERROR_out_of_memory) || (errcode == UMFPACK_ERROR_file_IO))
278  {
279  matrixIsLoaded_ = false;
280  setMatrix(mat_);
281  saveDecomposition(file);
282  }
283  else
284  {
285  matrixIsLoaded_ = true;
286  std::cout << "UMFPack decomposition successfully loaded from " << file << std::endl;
287  }
288  }
289 
296  UMFPack(const char* file, int verbose=0)
297  {
298  //check whether T is a supported type
299  static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
300  "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
301  Caller::defaults(UMF_Control);
302  int errcode = Caller::load_numeric(&UMF_Numeric, const_cast<char*>(file));
303  if (errcode == UMFPACK_ERROR_out_of_memory)
304  DUNE_THROW(Dune::Exception, "ran out of memory while loading UMFPack decomposition");
305  if (errcode == UMFPACK_ERROR_file_IO)
306  DUNE_THROW(Dune::Exception, "IO error while loading UMFPack decomposition");
307  matrixIsLoaded_ = true;
308  std::cout << "UMFPack decomposition successfully loaded from " << file << std::endl;
309  setVerbosity(verbose);
310  }
311 
312  virtual ~UMFPack()
313  {
314  if ((umfpackMatrix_.N() + umfpackMatrix_.M() > 0) || matrixIsLoaded_)
315  free();
316  }
317 
322  {
323  if (umfpackMatrix_.N() != b.dim())
324  DUNE_THROW(Dune::ISTLError, "Size of right-hand-side vector b does not match the number of matrix rows!");
325  if (umfpackMatrix_.M() != x.dim())
326  DUNE_THROW(Dune::ISTLError, "Size of solution vector x does not match the number of matrix columns!");
327 
328  double UMF_Apply_Info[UMFPACK_INFO];
329  Caller::solve(UMFPACK_A,
330  umfpackMatrix_.getColStart(),
331  umfpackMatrix_.getRowIndex(),
332  umfpackMatrix_.getValues(),
333  reinterpret_cast<double*>(&x[0]),
334  reinterpret_cast<double*>(&b[0]),
335  UMF_Numeric,
336  UMF_Control,
337  UMF_Apply_Info);
338 
339  //this is a direct solver
340  res.iterations = 1;
341  res.converged = true;
342  res.elapsed = UMF_Apply_Info[UMFPACK_SOLVE_WALLTIME];
343 
344  printOnApply(UMF_Apply_Info);
345  }
346 
350  virtual void apply (domain_type& x, range_type& b, double reduction, InverseOperatorResult& res)
351  {
352  DUNE_UNUSED_PARAMETER(reduction);
353  apply(x,b,res);
354  }
355 
361  void apply(T* x, T* b)
362  {
363  double UMF_Apply_Info[UMFPACK_INFO];
364  Caller::solve(UMFPACK_A,
365  umfpackMatrix_.getColStart(),
366  umfpackMatrix_.getRowIndex(),
367  umfpackMatrix_.getValues(),
368  x,
369  b,
370  UMF_Numeric,
371  UMF_Control,
372  UMF_Apply_Info);
373  printOnApply(UMF_Apply_Info);
374  }
375 
387  void setOption(unsigned int option, double value)
388  {
389  if (option >= UMFPACK_CONTROL)
390  DUNE_THROW(RangeError, "Requested non-existing UMFPack option");
391 
392  UMF_Control[option] = value;
393  }
394 
398  void saveDecomposition(const char* file)
399  {
400  int errcode = Caller::save_numeric(UMF_Numeric, const_cast<char*>(file));
401  if (errcode != UMFPACK_OK)
402  DUNE_THROW(Dune::Exception,"IO ERROR while trying to save UMFPack decomposition");
403  }
404 
406  void setMatrix(const Matrix& matrix)
407  {
408  if ((umfpackMatrix_.N() + umfpackMatrix_.M() > 0) || matrixIsLoaded_)
409  free();
410  umfpackMatrix_ = matrix;
411  decompose();
412  }
413 
414  template<class S>
415  void setSubMatrix(const Matrix& _mat, const S& rowIndexSet)
416  {
417  if ((umfpackMatrix_.N() + umfpackMatrix_.M() > 0) || matrixIsLoaded_)
418  free();
419  umfpackMatrix_.setMatrix(_mat,rowIndexSet);
420  decompose();
421  }
422 
430  void setVerbosity(int v)
431  {
432  verbosity_ = v;
433  // set the verbosity level in UMFPack
434  if (verbosity_ == 0)
435  UMF_Control[UMFPACK_PRL] = 1;
436  if (verbosity_ == 1)
437  UMF_Control[UMFPACK_PRL] = 2;
438  if (verbosity_ == 2)
439  UMF_Control[UMFPACK_PRL] = 4;
440  }
441 
447  {
448  return UMF_Numeric;
449  }
450 
456  {
457  return umfpackMatrix_;
458  }
459 
464  void free()
465  {
466  if (!matrixIsLoaded_)
467  {
468  Caller::free_symbolic(&UMF_Symbolic);
469  umfpackMatrix_.free();
470  }
471  Caller::free_numeric(&UMF_Numeric);
472  matrixIsLoaded_ = false;
473  }
474 
475  const char* name() { return "UMFPACK"; }
476 
477  private:
478  typedef typename Dune::UMFPackMethodChooser<T> Caller;
479 
480  template<class M,class X, class TM, class TD, class T1>
481  friend class SeqOverlappingSchwarz;
482  friend struct SeqOverlappingSchwarzAssemblerHelper<UMFPack<Matrix>,true>;
483 
485  void decompose()
486  {
487  double UMF_Decomposition_Info[UMFPACK_INFO];
488  Caller::symbolic(static_cast<int>(umfpackMatrix_.N()),
489  static_cast<int>(umfpackMatrix_.N()),
490  umfpackMatrix_.getColStart(),
491  umfpackMatrix_.getRowIndex(),
492  reinterpret_cast<double*>(umfpackMatrix_.getValues()),
493  &UMF_Symbolic,
494  UMF_Control,
495  UMF_Decomposition_Info);
496  Caller::numeric(umfpackMatrix_.getColStart(),
497  umfpackMatrix_.getRowIndex(),
498  reinterpret_cast<double*>(umfpackMatrix_.getValues()),
499  UMF_Symbolic,
500  &UMF_Numeric,
501  UMF_Control,
502  UMF_Decomposition_Info);
503  Caller::report_status(UMF_Control,UMF_Decomposition_Info[UMFPACK_STATUS]);
504  if (verbosity_ == 1)
505  {
506  std::cout << "[UMFPack Decomposition]" << std::endl;
507  std::cout << "Wallclock Time taken: " << UMF_Decomposition_Info[UMFPACK_NUMERIC_WALLTIME] << " (CPU Time: " << UMF_Decomposition_Info[UMFPACK_NUMERIC_TIME] << ")" << std::endl;
508  std::cout << "Flops taken: " << UMF_Decomposition_Info[UMFPACK_FLOPS] << std::endl;
509  std::cout << "Peak Memory Usage: " << UMF_Decomposition_Info[UMFPACK_PEAK_MEMORY]*UMF_Decomposition_Info[UMFPACK_SIZE_OF_UNIT] << " bytes" << std::endl;
510  std::cout << "Condition number estimate: " << 1./UMF_Decomposition_Info[UMFPACK_RCOND] << std::endl;
511  std::cout << "Numbers of non-zeroes in decomposition: L: " << UMF_Decomposition_Info[UMFPACK_LNZ] << " U: " << UMF_Decomposition_Info[UMFPACK_UNZ] << std::endl;
512  }
513  if (verbosity_ == 2)
514  {
515  Caller::report_info(UMF_Control,UMF_Decomposition_Info);
516  }
517  }
518 
519  void printOnApply(double* UMF_Info)
520  {
521  Caller::report_status(UMF_Control,UMF_Info[UMFPACK_STATUS]);
522  if (verbosity_ > 0)
523  {
524  std::cout << "[UMFPack Solve]" << std::endl;
525  std::cout << "Wallclock Time: " << UMF_Info[UMFPACK_SOLVE_WALLTIME] << " (CPU Time: " << UMF_Info[UMFPACK_SOLVE_TIME] << ")" << std::endl;
526  std::cout << "Flops Taken: " << UMF_Info[UMFPACK_SOLVE_FLOPS] << std::endl;
527  std::cout << "Iterative Refinement steps taken: " << UMF_Info[UMFPACK_IR_TAKEN] << std::endl;
528  std::cout << "Error Estimate: " << UMF_Info[UMFPACK_OMEGA1] << " resp. " << UMF_Info[UMFPACK_OMEGA2] << std::endl;
529  }
530  }
531 
532  UMFPackMatrix umfpackMatrix_;
533  bool matrixIsLoaded_;
534  int verbosity_;
535  void *UMF_Symbolic;
536  void *UMF_Numeric;
537  double UMF_Control[UMFPACK_CONTROL];
538  };
539 
540  template<typename T, typename A, int n, int m>
541  struct IsDirectSolver<UMFPack<BCRSMatrix<FieldMatrix<T,n,m>,A> > >
542  {
543  enum { value=true};
544  };
545 
546  template<typename T, typename A, int n, int m>
547  struct StoresColumnCompressed<UMFPack<BCRSMatrix<FieldMatrix<T,n,m>,A> > >
548  {
549  enum { value = true };
550  };
551 }
552 
553 #endif // HAVE_SUITESPARSE_UMFPACK
554 
555 #endif //DUNE_ISTL_UMFPACK_HH
Implementation of the BCRSMatrix class.
A sparse block matrix with compressed row storage.
Definition: bcrsmatrix.hh:423
A vector of blocks with memory management.
Definition: bvector.hh:313
Base class for Dune-Exceptions.
Definition: exceptions.hh:94
A dense n x m matrix.
Definition: fmatrix.hh:68
vector space out of a tensor product of fields.
Definition: fvector.hh:93
derive error class from the base class in common
Definition: istlexception.hh:16
Abstract base class for all solvers.
Definition: solver.hh:79
A generic dynamic dense matrix.
Definition: matrix.hh:555
Default exception class for range errors.
Definition: exceptions.hh:252
Use the UMFPack package to directly solve linear systems – empty default class.
Definition: umfpack.hh:50
size_type dim() const
dimension of the vector space
Definition: bvector.hh:283
A few common exception classes.
Implements a matrix constructed from a given type representing a field and compile-time given number ...
Implements a vector constructed from a given type representing a field and a compile-time given size.
#define DUNE_THROW(E, m)
Definition: exceptions.hh:216
UMFPack(const Matrix &matrix, int verbose, bool)
Constructor for compatibility with SuperLU standard constructor.
Definition: umfpack.hh:239
void * getFactorization()
Return the matrix factorization.
Definition: umfpack.hh:446
void setVerbosity(int v)
sets the verbosity level for the UMFPack solver
Definition: umfpack.hh:430
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:269
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:207
UMFPackMatrix & getInternalMatrix()
Return the column compress matrix from UMFPack.
Definition: umfpack.hh:455
void setMatrix(const Matrix &matrix)
Initialize data from given matrix.
Definition: umfpack.hh:406
virtual void apply(domain_type &x, range_type &b, InverseOperatorResult &res)
Apply inverse operator,.
Definition: umfpack.hh:321
void free()
free allocated space.
Definition: umfpack.hh:464
ColCompMatrixInitializer< BCRSMatrix< FieldMatrix< T, n, m >, A > > MatrixInitializer
Type of an associated initializer class.
Definition: umfpack.hh:203
UMFPack(const Matrix &matrix, int verbose=0)
Construct a solver object from a BCRSMatrix.
Definition: umfpack.hh:221
void setOption(unsigned int option, double value)
Set UMFPack-specific options.
Definition: umfpack.hh:387
Dune::ColCompMatrix< Matrix > UMFPackMatrix
The corresponding SuperLU Matrix type.
Definition: umfpack.hh:201
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:211
void apply(T *x, T *b)
additional apply method with c-arrays in analogy to superlu
Definition: umfpack.hh:361
UMFPack()
default constructor
Definition: umfpack.hh:251
UMFPack(const char *file, int verbose=0)
try loading a decomposition from file
Definition: umfpack.hh:296
virtual void apply(domain_type &x, range_type &b, double reduction, InverseOperatorResult &res)
apply inverse operator, with given convergence criteria.
Definition: umfpack.hh:350
void saveDecomposition(const char *file)
saves a decomposition to a file
Definition: umfpack.hh:398
Dune namespace.
Definition: alignment.hh:11
Implementations of the inverse operator interface.
Templates characterizing the type of a solver.
Inititializer for the ColCompMatrix as needed by OverlappingSchwarz.
Definition: colcompmatrix.hh:154
Statistics about the application of an inverse operator.
Definition: solver.hh:32
double elapsed
Elapsed time in seconds.
Definition: solver.hh:62
int iterations
Number of iterations.
Definition: solver.hh:50
bool converged
True if convergence criterion has been met.
Definition: solver.hh:56
Definition of the DUNE_UNUSED macro for the case that config.h is not available.
#define DUNE_UNUSED_PARAMETER(parm)
A macro to mark intentionally unused function parameters with.
Definition: unused.hh:18
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (May 9, 22:29, 2024)