DUNE PDELab (2.8)

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