Dune Core Modules (2.9.0)

umfpack.hh
Go to the documentation of this file.
1 // SPDX-FileCopyrightText: Copyright (C) DUNE Project contributors, see file LICENSE.md in module root
2 // SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
3 // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
4 // vi: set et ts=4 sw=2 sts=2:
5 #ifndef DUNE_ISTL_UMFPACK_HH
6 #define DUNE_ISTL_UMFPACK_HH
7 
8 #if HAVE_SUITESPARSE_UMFPACK || defined DOXYGEN
9 
10 #include<complex>
11 #include<type_traits>
12 
13 #include<umfpack.h>
14 
16 #include<dune/common/fmatrix.hh>
17 #include<dune/common/fvector.hh>
18 #include<dune/istl/bccsmatrixinitializer.hh>
20 #include<dune/istl/solvers.hh>
22 #include <dune/istl/solverfactory.hh>
23 
24 
25 
26 namespace Dune {
38  // FORWARD DECLARATIONS
39  template<class M, class T, class TM, class TD, class TA>
40  class SeqOverlappingSchwarz;
41 
42  template<class T, bool tag>
43  struct SeqOverlappingSchwarzAssemblerHelper;
44 
45  // wrapper class for C-Function Calls in the backend. Choose the right function namespace
46  // depending on the template parameter used.
47  template<typename T>
48  struct UMFPackMethodChooser
49  {
50  static constexpr bool valid = false ;
51  };
52 
53  template<>
54  struct UMFPackMethodChooser<double>
55  {
56  static constexpr bool valid = true ;
57 
58  template<typename... A>
59  static void defaults(A... args)
60  {
61  umfpack_dl_defaults(args...);
62  }
63  template<typename... A>
64  static void free_numeric(A... args)
65  {
66  umfpack_dl_free_numeric(args...);
67  }
68  template<typename... A>
69  static void free_symbolic(A... args)
70  {
71  umfpack_dl_free_symbolic(args...);
72  }
73  template<typename... A>
74  static int load_numeric(A... args)
75  {
76  return umfpack_dl_load_numeric(args...);
77  }
78  template<typename... A>
79  static void numeric(A... args)
80  {
81  umfpack_dl_numeric(args...);
82  }
83  template<typename... A>
84  static void report_info(A... args)
85  {
86  umfpack_dl_report_info(args...);
87  }
88  template<typename... A>
89  static void report_status(A... args)
90  {
91  umfpack_dl_report_status(args...);
92  }
93  template<typename... A>
94  static int save_numeric(A... args)
95  {
96  return umfpack_dl_save_numeric(args...);
97  }
98  template<typename... A>
99  static void solve(A... args)
100  {
101  umfpack_dl_solve(args...);
102  }
103  template<typename... A>
104  static void symbolic(A... args)
105  {
106  umfpack_dl_symbolic(args...);
107  }
108  };
109 
110  template<>
111  struct UMFPackMethodChooser<std::complex<double> >
112  {
113  static constexpr bool valid = true ;
114 
115  template<typename... A>
116  static void defaults(A... args)
117  {
118  umfpack_zl_defaults(args...);
119  }
120  template<typename... A>
121  static void free_numeric(A... args)
122  {
123  umfpack_zl_free_numeric(args...);
124  }
125  template<typename... A>
126  static void free_symbolic(A... args)
127  {
128  umfpack_zl_free_symbolic(args...);
129  }
130  template<typename... A>
131  static int load_numeric(A... args)
132  {
133  return umfpack_zl_load_numeric(args...);
134  }
135  template<typename... A>
136  static void numeric(const long int* cs, const long int* ri, const double* val, A... args)
137  {
138  umfpack_zl_numeric(cs,ri,val,NULL,args...);
139  }
140  template<typename... A>
141  static void report_info(A... args)
142  {
143  umfpack_zl_report_info(args...);
144  }
145  template<typename... A>
146  static void report_status(A... args)
147  {
148  umfpack_zl_report_status(args...);
149  }
150  template<typename... A>
151  static int save_numeric(A... args)
152  {
153  return umfpack_zl_save_numeric(args...);
154  }
155  template<typename... A>
156  static void solve(long int m, const long int* cs, const long int* ri, std::complex<double>* val, double* x, const double* b,A... args)
157  {
158  const double* cval = reinterpret_cast<const double*>(val);
159  umfpack_zl_solve(m,cs,ri,cval,NULL,x,NULL,b,NULL,args...);
160  }
161  template<typename... A>
162  static void symbolic(long int m, long int n, const long int* cs, const long int* ri, const double* val, A... args)
163  {
164  umfpack_zl_symbolic(m,n,cs,ri,val,NULL,args...);
165  }
166  };
167 
168  namespace Impl
169  {
170  template<class M>
171  struct UMFPackVectorChooser
172  {};
173 
174  template<typename T, typename A, int n, int m>
175  struct UMFPackVectorChooser<BCRSMatrix<FieldMatrix<T,n,m>,A > >
176  {
178  using domain_type = BlockVector<
179  FieldVector<T,m>,
180  typename std::allocator_traits<A>::template rebind_alloc<FieldVector<T,m> > >;
182  using range_type = BlockVector<
183  FieldVector<T,n>,
184  typename std::allocator_traits<A>::template rebind_alloc<FieldVector<T,n> > >;
185  };
186 
187  template<typename T, typename A>
188  struct UMFPackVectorChooser<BCRSMatrix<T,A> >
189  {
191  using domain_type = BlockVector<T, A>;
193  using range_type = BlockVector<T, A>;
194  };
195  }
196 
210  template<typename M>
211  class UMFPack
212  : public InverseOperator<
213  typename Impl::UMFPackVectorChooser<M>::domain_type,
214  typename Impl::UMFPackVectorChooser<M>::range_type >
215  {
216  using T = typename M::field_type;
217 
218  public:
220  using Matrix = M;
221  using matrix_type = M;
223  typedef ISTL::Impl::BCCSMatrix<typename Matrix::field_type, long int> UMFPackMatrix;
225  typedef ISTL::Impl::BCCSMatrixInitializer<M, long int> MatrixInitializer;
227  using domain_type = typename Impl::UMFPackVectorChooser<M>::domain_type;
229  using range_type = typename Impl::UMFPackVectorChooser<M>::range_type;
230 
233  {
234  return SolverCategory::Category::sequential;
235  }
236 
245  UMFPack(const Matrix& matrix, int verbose=0) : matrixIsLoaded_(false)
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  setVerbosity(verbose);
252  setMatrix(matrix);
253  }
254 
263  UMFPack(const Matrix& matrix, int verbose, bool) : matrixIsLoaded_(false)
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  setMatrix(matrix);
271  }
272 
282  UMFPack(const Matrix& mat_, const ParameterTree& config)
283  : UMFPack(mat_, config.get<int>("verbose", 0))
284  {}
285 
288  UMFPack() : matrixIsLoaded_(false), verbosity_(0)
289  {
290  //check whether T is a supported type
291  static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
292  "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
293  Caller::defaults(UMF_Control);
294  }
295 
306  UMFPack(const Matrix& mat_, const char* file, int verbose=0)
307  {
308  //check whether T is a supported type
309  static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
310  "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
311  Caller::defaults(UMF_Control);
312  setVerbosity(verbose);
313  int errcode = Caller::load_numeric(&UMF_Numeric, const_cast<char*>(file));
314  if ((errcode == UMFPACK_ERROR_out_of_memory) || (errcode == UMFPACK_ERROR_file_IO))
315  {
316  matrixIsLoaded_ = false;
317  setMatrix(mat_);
318  saveDecomposition(file);
319  }
320  else
321  {
322  matrixIsLoaded_ = true;
323  std::cout << "UMFPack decomposition successfully loaded from " << file << std::endl;
324  }
325  }
326 
333  UMFPack(const char* file, int verbose=0)
334  {
335  //check whether T is a supported type
336  static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
337  "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
338  Caller::defaults(UMF_Control);
339  int errcode = Caller::load_numeric(&UMF_Numeric, const_cast<char*>(file));
340  if (errcode == UMFPACK_ERROR_out_of_memory)
341  DUNE_THROW(Dune::Exception, "ran out of memory while loading UMFPack decomposition");
342  if (errcode == UMFPACK_ERROR_file_IO)
343  DUNE_THROW(Dune::Exception, "IO error while loading UMFPack decomposition");
344  matrixIsLoaded_ = true;
345  std::cout << "UMFPack decomposition successfully loaded from " << file << std::endl;
346  setVerbosity(verbose);
347  }
348 
349  virtual ~UMFPack()
350  {
351  if ((umfpackMatrix_.N() + umfpackMatrix_.M() > 0) || matrixIsLoaded_)
352  free();
353  }
354 
359  {
360  if (umfpackMatrix_.N() != b.dim())
361  DUNE_THROW(Dune::ISTLError, "Size of right-hand-side vector b does not match the number of matrix rows!");
362  if (umfpackMatrix_.M() != x.dim())
363  DUNE_THROW(Dune::ISTLError, "Size of solution vector x does not match the number of matrix columns!");
364  if (b.size() == 0)
365  return;
366 
367  double UMF_Apply_Info[UMFPACK_INFO];
368  Caller::solve(UMFPACK_A,
369  umfpackMatrix_.getColStart(),
370  umfpackMatrix_.getRowIndex(),
371  umfpackMatrix_.getValues(),
372  reinterpret_cast<double*>(&x[0]),
373  reinterpret_cast<double*>(&b[0]),
374  UMF_Numeric,
375  UMF_Control,
376  UMF_Apply_Info);
377 
378  //this is a direct solver
379  res.iterations = 1;
380  res.converged = true;
381  res.elapsed = UMF_Apply_Info[UMFPACK_SOLVE_WALLTIME];
382 
383  printOnApply(UMF_Apply_Info);
384  }
385 
389  virtual void apply (domain_type& x, range_type& b, [[maybe_unused]] double reduction, InverseOperatorResult& res)
390  {
391  apply(x,b,res);
392  }
393 
399  void apply(T* x, T* b)
400  {
401  double UMF_Apply_Info[UMFPACK_INFO];
402  Caller::solve(UMFPACK_A,
403  umfpackMatrix_.getColStart(),
404  umfpackMatrix_.getRowIndex(),
405  umfpackMatrix_.getValues(),
406  x,
407  b,
408  UMF_Numeric,
409  UMF_Control,
410  UMF_Apply_Info);
411  printOnApply(UMF_Apply_Info);
412  }
413 
425  void setOption(unsigned int option, double value)
426  {
427  if (option >= UMFPACK_CONTROL)
428  DUNE_THROW(RangeError, "Requested non-existing UMFPack option");
429 
430  UMF_Control[option] = value;
431  }
432 
436  void saveDecomposition(const char* file)
437  {
438  int errcode = Caller::save_numeric(UMF_Numeric, const_cast<char*>(file));
439  if (errcode != UMFPACK_OK)
440  DUNE_THROW(Dune::Exception,"IO ERROR while trying to save UMFPack decomposition");
441  }
442 
444  void setMatrix(const Matrix& matrix)
445  {
446  if ((umfpackMatrix_.N() + umfpackMatrix_.M() > 0) || matrixIsLoaded_)
447  free();
448  if (matrix.N() == 0 or matrix.M() == 0)
449  return;
450 
451  if (umfpackMatrix_.N() + umfpackMatrix_.M() + umfpackMatrix_.nonzeroes() != 0)
452  umfpackMatrix_.free();
453  umfpackMatrix_.setSize(MatrixDimension<Matrix>::rowdim(matrix),
454  MatrixDimension<Matrix>::coldim(matrix));
455  ISTL::Impl::BCCSMatrixInitializer<Matrix, long int> initializer(umfpackMatrix_);
456 
457  copyToBCCSMatrix(initializer, matrix);
458 
459  decompose();
460  }
461 
462  template<class S>
463  void setSubMatrix(const Matrix& _mat, const S& rowIndexSet)
464  {
465  if ((umfpackMatrix_.N() + umfpackMatrix_.M() > 0) || matrixIsLoaded_)
466  free();
467 
468  if (umfpackMatrix_.N() + umfpackMatrix_.M() + umfpackMatrix_.nonzeroes() != 0)
469  umfpackMatrix_.free();
470 
471  umfpackMatrix_.setSize(rowIndexSet.size()*MatrixDimension<Matrix>::rowdim(_mat) / _mat.N(),
472  rowIndexSet.size()*MatrixDimension<Matrix>::coldim(_mat) / _mat.M());
473  ISTL::Impl::BCCSMatrixInitializer<Matrix, long int> initializer(umfpackMatrix_);
474 
475  copyToBCCSMatrix(initializer, ISTL::Impl::MatrixRowSubset<Matrix,std::set<std::size_t> >(_mat,rowIndexSet));
476 
477  decompose();
478  }
479 
487  void setVerbosity(int v)
488  {
489  verbosity_ = v;
490  // set the verbosity level in UMFPack
491  if (verbosity_ == 0)
492  UMF_Control[UMFPACK_PRL] = 1;
493  if (verbosity_ == 1)
494  UMF_Control[UMFPACK_PRL] = 2;
495  if (verbosity_ == 2)
496  UMF_Control[UMFPACK_PRL] = 4;
497  }
498 
504  {
505  return UMF_Numeric;
506  }
507 
513  {
514  return umfpackMatrix_;
515  }
516 
521  void free()
522  {
523  if (!matrixIsLoaded_)
524  {
525  Caller::free_symbolic(&UMF_Symbolic);
526  umfpackMatrix_.free();
527  }
528  Caller::free_numeric(&UMF_Numeric);
529  matrixIsLoaded_ = false;
530  }
531 
532  const char* name() { return "UMFPACK"; }
533 
534  private:
535  typedef typename Dune::UMFPackMethodChooser<T> Caller;
536 
537  template<class Mat,class X, class TM, class TD, class T1>
538  friend class SeqOverlappingSchwarz;
539  friend struct SeqOverlappingSchwarzAssemblerHelper<UMFPack<Matrix>,true>;
540 
542  void decompose()
543  {
544  double UMF_Decomposition_Info[UMFPACK_INFO];
545  Caller::symbolic(static_cast<int>(umfpackMatrix_.N()),
546  static_cast<int>(umfpackMatrix_.N()),
547  umfpackMatrix_.getColStart(),
548  umfpackMatrix_.getRowIndex(),
549  reinterpret_cast<double*>(umfpackMatrix_.getValues()),
550  &UMF_Symbolic,
551  UMF_Control,
552  UMF_Decomposition_Info);
553  Caller::numeric(umfpackMatrix_.getColStart(),
554  umfpackMatrix_.getRowIndex(),
555  reinterpret_cast<double*>(umfpackMatrix_.getValues()),
556  UMF_Symbolic,
557  &UMF_Numeric,
558  UMF_Control,
559  UMF_Decomposition_Info);
560  Caller::report_status(UMF_Control,UMF_Decomposition_Info[UMFPACK_STATUS]);
561  if (verbosity_ == 1)
562  {
563  std::cout << "[UMFPack Decomposition]" << std::endl;
564  std::cout << "Wallclock Time taken: " << UMF_Decomposition_Info[UMFPACK_NUMERIC_WALLTIME] << " (CPU Time: " << UMF_Decomposition_Info[UMFPACK_NUMERIC_TIME] << ")" << std::endl;
565  std::cout << "Flops taken: " << UMF_Decomposition_Info[UMFPACK_FLOPS] << std::endl;
566  std::cout << "Peak Memory Usage: " << UMF_Decomposition_Info[UMFPACK_PEAK_MEMORY]*UMF_Decomposition_Info[UMFPACK_SIZE_OF_UNIT] << " bytes" << std::endl;
567  std::cout << "Condition number estimate: " << 1./UMF_Decomposition_Info[UMFPACK_RCOND] << std::endl;
568  std::cout << "Numbers of non-zeroes in decomposition: L: " << UMF_Decomposition_Info[UMFPACK_LNZ] << " U: " << UMF_Decomposition_Info[UMFPACK_UNZ] << std::endl;
569  }
570  if (verbosity_ == 2)
571  {
572  Caller::report_info(UMF_Control,UMF_Decomposition_Info);
573  }
574  }
575 
576  void printOnApply(double* UMF_Info)
577  {
578  Caller::report_status(UMF_Control,UMF_Info[UMFPACK_STATUS]);
579  if (verbosity_ > 0)
580  {
581  std::cout << "[UMFPack Solve]" << std::endl;
582  std::cout << "Wallclock Time: " << UMF_Info[UMFPACK_SOLVE_WALLTIME] << " (CPU Time: " << UMF_Info[UMFPACK_SOLVE_TIME] << ")" << std::endl;
583  std::cout << "Flops Taken: " << UMF_Info[UMFPACK_SOLVE_FLOPS] << std::endl;
584  std::cout << "Iterative Refinement steps taken: " << UMF_Info[UMFPACK_IR_TAKEN] << std::endl;
585  std::cout << "Error Estimate: " << UMF_Info[UMFPACK_OMEGA1] << " resp. " << UMF_Info[UMFPACK_OMEGA2] << std::endl;
586  }
587  }
588 
589  UMFPackMatrix umfpackMatrix_;
590  bool matrixIsLoaded_;
591  int verbosity_;
592  void *UMF_Symbolic;
593  void *UMF_Numeric;
594  double UMF_Control[UMFPACK_CONTROL];
595  };
596 
597  template<typename T, typename A, int n, int m>
598  struct IsDirectSolver<UMFPack<BCRSMatrix<FieldMatrix<T,n,m>,A> > >
599  {
600  enum { value=true};
601  };
602 
603  template<typename T, typename A>
604  struct StoresColumnCompressed<UMFPack<BCRSMatrix<T,A> > >
605  {
606  enum { value = true };
607  };
608 
609  struct UMFPackCreator {
610  template<class F,class=void> struct isValidBlock : std::false_type{};
611  template<class B> struct isValidBlock<B, std::enable_if_t<std::is_same<typename FieldTraits<B>::real_type,double>::value>> : std::true_type {};
612 
613  template<typename TL, typename M>
614  std::shared_ptr<Dune::InverseOperator<typename Dune::TypeListElement<1, TL>::type,
615  typename Dune::TypeListElement<2, TL>::type>>
616  operator() (TL /*tl*/, const M& mat, const Dune::ParameterTree& config,
617  std::enable_if_t<
618  isValidBlock<typename Dune::TypeListElement<1, TL>::type::block_type>::value,int> = 0) const
619  {
620  int verbose = config.get("verbose", 0);
621  return std::make_shared<Dune::UMFPack<M>>(mat,verbose);
622  }
623 
624  // second version with SFINAE to validate the template parameters of UMFPack
625  template<typename TL, typename M>
626  std::shared_ptr<Dune::InverseOperator<typename Dune::TypeListElement<1, TL>::type,
627  typename Dune::TypeListElement<2, TL>::type>>
628  operator() (TL /*tl*/, const M& /*mat*/, const Dune::ParameterTree& /*config*/,
629  std::enable_if_t<
630  !isValidBlock<typename Dune::TypeListElement<1, TL>::type::block_type>::value,int> = 0) const
631  {
632  DUNE_THROW(UnsupportedType,
633  "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
634  }
635  };
636  DUNE_REGISTER_DIRECT_SOLVER("umfpack",Dune::UMFPackCreator());
637 } // end namespace Dune
638 
639 #endif // HAVE_SUITESPARSE_UMFPACK
640 
641 #endif //DUNE_ISTL_UMFPACK_HH
Implementation of the BCRSMatrix class.
Base class for Dune-Exceptions.
Definition: exceptions.hh:96
derive error class from the base class in common
Definition: istlexception.hh:19
Abstract base class for all solvers.
Definition: solver.hh:99
A generic dynamic dense matrix.
Definition: matrix.hh:561
size_type M() const
Return the number of columns.
Definition: matrix.hh:700
size_type N() const
Return the number of rows.
Definition: matrix.hh:695
Hierarchical structure of string parameters.
Definition: parametertree.hh:37
std::string get(const std::string &key, const std::string &defaultValue) const
get value as string
Definition: parametertree.cc:185
Default exception class for range errors.
Definition: exceptions.hh:254
The UMFPack direct sparse solver.
Definition: umfpack.hh:215
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:218
void free()
free allocated space.
Definition: umfpack.hh:521
virtual void apply(domain_type &x, range_type &b, [[maybe_unused]] double reduction, InverseOperatorResult &res)
apply inverse operator, with given convergence criteria.
Definition: umfpack.hh:389
virtual void apply(domain_type &x, range_type &b, InverseOperatorResult &res)
Apply inverse operator,.
Definition: umfpack.hh:358
void * getFactorization()
Return the matrix factorization.
Definition: umfpack.hh:503
virtual SolverCategory::Category category() const
Category of the solver (see SolverCategory::Category)
Definition: umfpack.hh:232
void setMatrix(const Matrix &matrix)
Initialize data from given matrix.
Definition: umfpack.hh:444
UMFPack(const Matrix &mat_, const ParameterTree &config)
Construct a solver object from a matrix.
Definition: umfpack.hh:282
ISTL::Impl::BCCSMatrix< typename Matrix::field_type, long int > UMFPackMatrix
The corresponding UMFPack matrix type.
Definition: umfpack.hh:223
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:306
typename Impl::UMFPackVectorChooser< M >::range_type range_type
The type of the range of the solver.
Definition: umfpack.hh:229
UMFPack()
default constructor
Definition: umfpack.hh:288
UMFPackMatrix & getInternalMatrix()
Return the column compress matrix from UMFPack.
Definition: umfpack.hh:512
ISTL::Impl::BCCSMatrixInitializer< M, long int > MatrixInitializer
Type of an associated initializer class.
Definition: umfpack.hh:225
void apply(T *x, T *b)
additional apply method with c-arrays in analogy to superlu
Definition: umfpack.hh:399
void setVerbosity(int v)
sets the verbosity level for the UMFPack solver
Definition: umfpack.hh:487
UMFPack(const char *file, int verbose=0)
try loading a decomposition from file
Definition: umfpack.hh:333
void saveDecomposition(const char *file)
saves a decomposition to a file
Definition: umfpack.hh:436
typename Impl::UMFPackVectorChooser< M >::domain_type domain_type
The type of the domain of the solver.
Definition: umfpack.hh:227
UMFPack(const Matrix &matrix, int verbose=0)
Construct a solver object from a matrix.
Definition: umfpack.hh:245
void setOption(unsigned int option, double value)
Set UMFPack-specific options.
Definition: umfpack.hh:425
M Matrix
The matrix type.
Definition: umfpack.hh:220
UMFPack(const Matrix &matrix, int verbose, bool)
Constructor for compatibility with SuperLU standard constructor.
Definition: umfpack.hh:263
Dune namespace.
Definition: alignedallocator.hh:13
Implementations of the inverse operator interface.
Templates characterizing the type of a solver.
Statistics about the application of an inverse operator.
Definition: solver.hh:48
double elapsed
Elapsed time in seconds.
Definition: solver.hh:82
int iterations
Number of iterations.
Definition: solver.hh:67
bool converged
True if convergence criterion has been met.
Definition: solver.hh:73
Category
Definition: solvercategory.hh:23
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.80.0 (May 1, 22:29, 2024)