Dune Core Modules (2.7.1)

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
20#include <dune/istl/solverfactory.hh>
21
22#include"colcompmatrix.hh"
23
24
25namespace Dune {
37 // FORWARD DECLARATIONS
38 template<class M, class T, class TM, class TD, class TA>
39 class SeqOverlappingSchwarz;
40
41 template<class T, bool tag>
42 struct SeqOverlappingSchwarzAssemblerHelper;
43
44 // wrapper class for C-Function Calls in the backend. Choose the right function namespace
45 // depending on the template parameter used.
46 template<typename T>
47 struct UMFPackMethodChooser
48 {
49 static constexpr bool valid = false ;
50 };
51
52 template<>
53 struct UMFPackMethodChooser<double>
54 {
55 static constexpr bool valid = true ;
56
57 template<typename... A>
58 static void defaults(A... args)
59 {
60 umfpack_dl_defaults(args...);
61 }
62 template<typename... A>
63 static void free_numeric(A... args)
64 {
65 umfpack_dl_free_numeric(args...);
66 }
67 template<typename... A>
68 static void free_symbolic(A... args)
69 {
70 umfpack_dl_free_symbolic(args...);
71 }
72 template<typename... A>
73 static int load_numeric(A... args)
74 {
75 return umfpack_dl_load_numeric(args...);
76 }
77 template<typename... A>
78 static void numeric(A... args)
79 {
80 umfpack_dl_numeric(args...);
81 }
82 template<typename... A>
83 static void report_info(A... args)
84 {
85 umfpack_dl_report_info(args...);
86 }
87 template<typename... A>
88 static void report_status(A... args)
89 {
90 umfpack_dl_report_status(args...);
91 }
92 template<typename... A>
93 static int save_numeric(A... args)
94 {
95 return umfpack_dl_save_numeric(args...);
96 }
97 template<typename... A>
98 static void solve(A... args)
99 {
100 umfpack_dl_solve(args...);
101 }
102 template<typename... A>
103 static void symbolic(A... args)
104 {
105 umfpack_dl_symbolic(args...);
106 }
107 };
108
109 template<>
110 struct UMFPackMethodChooser<std::complex<double> >
111 {
112 static constexpr bool valid = true ;
113
114 template<typename... A>
115 static void defaults(A... args)
116 {
117 umfpack_zl_defaults(args...);
118 }
119 template<typename... A>
120 static void free_numeric(A... args)
121 {
122 umfpack_zl_free_numeric(args...);
123 }
124 template<typename... A>
125 static void free_symbolic(A... args)
126 {
127 umfpack_zl_free_symbolic(args...);
128 }
129 template<typename... A>
130 static int load_numeric(A... args)
131 {
132 return umfpack_zl_load_numeric(args...);
133 }
134 template<typename... A>
135 static void numeric(const long int* cs, const long int* ri, const double* val, A... args)
136 {
137 umfpack_zl_numeric(cs,ri,val,NULL,args...);
138 }
139 template<typename... A>
140 static void report_info(A... args)
141 {
142 umfpack_zl_report_info(args...);
143 }
144 template<typename... A>
145 static void report_status(A... args)
146 {
147 umfpack_zl_report_status(args...);
148 }
149 template<typename... A>
150 static int save_numeric(A... args)
151 {
152 return umfpack_zl_save_numeric(args...);
153 }
154 template<typename... A>
155 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)
156 {
157 const double* cval = reinterpret_cast<const double*>(val);
158 umfpack_zl_solve(m,cs,ri,cval,NULL,x,NULL,b,NULL,args...);
159 }
160 template<typename... A>
161 static void symbolic(long int m, long int n, const long int* cs, const long int* ri, const double* val, A... args)
162 {
163 umfpack_zl_symbolic(m,n,cs,ri,val,NULL,args...);
164 }
165 };
166
167 namespace Impl
168 {
169 template<class M>
170 struct UMFPackVectorChooser
171 {};
172
173 template<typename T, typename A, int n, int m>
174 struct UMFPackVectorChooser<BCRSMatrix<FieldMatrix<T,n,m>,A > >
175 {
177 using domain_type = BlockVector<
178 FieldVector<T,m>,
179 typename std::allocator_traits<A>::template rebind_alloc<FieldVector<T,m> > >;
181 using range_type = BlockVector<
182 FieldVector<T,n>,
183 typename std::allocator_traits<A>::template rebind_alloc<FieldVector<T,n> > >;
184 };
185
186 template<typename T, typename A>
187 struct UMFPackVectorChooser<BCRSMatrix<T,A> >
188 {
190 using domain_type = BlockVector<T, A>;
192 using range_type = BlockVector<T, A>;
193 };
194 }
195
209 template<typename M>
211 : public InverseOperator<
212 typename Impl::UMFPackVectorChooser<M>::domain_type,
213 typename Impl::UMFPackVectorChooser<M>::range_type >
214 {
215 using T = typename M::field_type;
216
217 public:
219 using Matrix = M;
220 using matrix_type = M;
226 using domain_type = typename Impl::UMFPackVectorChooser<M>::domain_type;
228 using range_type = typename Impl::UMFPackVectorChooser<M>::range_type;
229
232 {
233 return SolverCategory::Category::sequential;
234 }
235
244 UMFPack(const Matrix& matrix, int verbose=0) : matrixIsLoaded_(false)
245 {
246 //check whether T is a supported type
247 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
248 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
249 Caller::defaults(UMF_Control);
250 setVerbosity(verbose);
251 setMatrix(matrix);
252 }
253
262 UMFPack(const Matrix& matrix, int verbose, bool) : matrixIsLoaded_(false)
263 {
264 //check whether T is a supported type
265 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
266 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
267 Caller::defaults(UMF_Control);
268 setVerbosity(verbose);
269 setMatrix(matrix);
270 }
271
274 UMFPack() : matrixIsLoaded_(false), verbosity_(0)
275 {
276 //check whether T is a supported type
277 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
278 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
279 Caller::defaults(UMF_Control);
280 }
281
292 UMFPack(const Matrix& mat_, const char* file, int verbose=0)
293 {
294 //check whether T is a supported type
295 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
296 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
297 Caller::defaults(UMF_Control);
298 setVerbosity(verbose);
299 int errcode = Caller::load_numeric(&UMF_Numeric, const_cast<char*>(file));
300 if ((errcode == UMFPACK_ERROR_out_of_memory) || (errcode == UMFPACK_ERROR_file_IO))
301 {
302 matrixIsLoaded_ = false;
303 setMatrix(mat_);
304 saveDecomposition(file);
305 }
306 else
307 {
308 matrixIsLoaded_ = true;
309 std::cout << "UMFPack decomposition successfully loaded from " << file << std::endl;
310 }
311 }
312
319 UMFPack(const char* file, int verbose=0)
320 {
321 //check whether T is a supported type
322 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
323 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
324 Caller::defaults(UMF_Control);
325 int errcode = Caller::load_numeric(&UMF_Numeric, const_cast<char*>(file));
326 if (errcode == UMFPACK_ERROR_out_of_memory)
327 DUNE_THROW(Dune::Exception, "ran out of memory while loading UMFPack decomposition");
328 if (errcode == UMFPACK_ERROR_file_IO)
329 DUNE_THROW(Dune::Exception, "IO error while loading UMFPack decomposition");
330 matrixIsLoaded_ = true;
331 std::cout << "UMFPack decomposition successfully loaded from " << file << std::endl;
332 setVerbosity(verbose);
333 }
334
335 virtual ~UMFPack()
336 {
337 if ((umfpackMatrix_.N() + umfpackMatrix_.M() > 0) || matrixIsLoaded_)
338 free();
339 }
340
345 {
346 if (umfpackMatrix_.N() != b.dim())
347 DUNE_THROW(Dune::ISTLError, "Size of right-hand-side vector b does not match the number of matrix rows!");
348 if (umfpackMatrix_.M() != x.dim())
349 DUNE_THROW(Dune::ISTLError, "Size of solution vector x does not match the number of matrix columns!");
350 if (b.size() == 0)
351 return;
352
353 double UMF_Apply_Info[UMFPACK_INFO];
354 Caller::solve(UMFPACK_A,
355 umfpackMatrix_.getColStart(),
356 umfpackMatrix_.getRowIndex(),
357 umfpackMatrix_.getValues(),
358 reinterpret_cast<double*>(&x[0]),
359 reinterpret_cast<double*>(&b[0]),
360 UMF_Numeric,
361 UMF_Control,
362 UMF_Apply_Info);
363
364 //this is a direct solver
365 res.iterations = 1;
366 res.converged = true;
367 res.elapsed = UMF_Apply_Info[UMFPACK_SOLVE_WALLTIME];
368
369 printOnApply(UMF_Apply_Info);
370 }
371
375 virtual void apply (domain_type& x, range_type& b, double reduction, InverseOperatorResult& res)
376 {
377 DUNE_UNUSED_PARAMETER(reduction);
378 apply(x,b,res);
379 }
380
386 void apply(T* x, T* b)
387 {
388 double UMF_Apply_Info[UMFPACK_INFO];
389 Caller::solve(UMFPACK_A,
390 umfpackMatrix_.getColStart(),
391 umfpackMatrix_.getRowIndex(),
392 umfpackMatrix_.getValues(),
393 x,
394 b,
395 UMF_Numeric,
396 UMF_Control,
397 UMF_Apply_Info);
398 printOnApply(UMF_Apply_Info);
399 }
400
412 void setOption(unsigned int option, double value)
413 {
414 if (option >= UMFPACK_CONTROL)
415 DUNE_THROW(RangeError, "Requested non-existing UMFPack option");
416
417 UMF_Control[option] = value;
418 }
419
423 void saveDecomposition(const char* file)
424 {
425 int errcode = Caller::save_numeric(UMF_Numeric, const_cast<char*>(file));
426 if (errcode != UMFPACK_OK)
427 DUNE_THROW(Dune::Exception,"IO ERROR while trying to save UMFPack decomposition");
428 }
429
431 void setMatrix(const Matrix& matrix)
432 {
433 if ((umfpackMatrix_.N() + umfpackMatrix_.M() > 0) || matrixIsLoaded_)
434 free();
435 if (matrix.N() == 0 or matrix.M() == 0)
436 return;
437 umfpackMatrix_ = matrix;
438 decompose();
439 }
440
441 template<class S>
442 void setSubMatrix(const Matrix& _mat, const S& rowIndexSet)
443 {
444 if ((umfpackMatrix_.N() + umfpackMatrix_.M() > 0) || matrixIsLoaded_)
445 free();
446 umfpackMatrix_.setMatrix(_mat,rowIndexSet);
447 decompose();
448 }
449
457 void setVerbosity(int v)
458 {
459 verbosity_ = v;
460 // set the verbosity level in UMFPack
461 if (verbosity_ == 0)
462 UMF_Control[UMFPACK_PRL] = 1;
463 if (verbosity_ == 1)
464 UMF_Control[UMFPACK_PRL] = 2;
465 if (verbosity_ == 2)
466 UMF_Control[UMFPACK_PRL] = 4;
467 }
468
474 {
475 return UMF_Numeric;
476 }
477
483 {
484 return umfpackMatrix_;
485 }
486
491 void free()
492 {
493 if (!matrixIsLoaded_)
494 {
495 Caller::free_symbolic(&UMF_Symbolic);
496 umfpackMatrix_.free();
497 }
498 Caller::free_numeric(&UMF_Numeric);
499 matrixIsLoaded_ = false;
500 }
501
502 const char* name() { return "UMFPACK"; }
503
504 private:
505 typedef typename Dune::UMFPackMethodChooser<T> Caller;
506
507 template<class Mat,class X, class TM, class TD, class T1>
508 friend class SeqOverlappingSchwarz;
509 friend struct SeqOverlappingSchwarzAssemblerHelper<UMFPack<Matrix>,true>;
510
512 void decompose()
513 {
514 double UMF_Decomposition_Info[UMFPACK_INFO];
515 Caller::symbolic(static_cast<int>(umfpackMatrix_.N()),
516 static_cast<int>(umfpackMatrix_.N()),
517 umfpackMatrix_.getColStart(),
518 umfpackMatrix_.getRowIndex(),
519 reinterpret_cast<double*>(umfpackMatrix_.getValues()),
520 &UMF_Symbolic,
521 UMF_Control,
522 UMF_Decomposition_Info);
523 Caller::numeric(umfpackMatrix_.getColStart(),
524 umfpackMatrix_.getRowIndex(),
525 reinterpret_cast<double*>(umfpackMatrix_.getValues()),
526 UMF_Symbolic,
527 &UMF_Numeric,
528 UMF_Control,
529 UMF_Decomposition_Info);
530 Caller::report_status(UMF_Control,UMF_Decomposition_Info[UMFPACK_STATUS]);
531 if (verbosity_ == 1)
532 {
533 std::cout << "[UMFPack Decomposition]" << std::endl;
534 std::cout << "Wallclock Time taken: " << UMF_Decomposition_Info[UMFPACK_NUMERIC_WALLTIME] << " (CPU Time: " << UMF_Decomposition_Info[UMFPACK_NUMERIC_TIME] << ")" << std::endl;
535 std::cout << "Flops taken: " << UMF_Decomposition_Info[UMFPACK_FLOPS] << std::endl;
536 std::cout << "Peak Memory Usage: " << UMF_Decomposition_Info[UMFPACK_PEAK_MEMORY]*UMF_Decomposition_Info[UMFPACK_SIZE_OF_UNIT] << " bytes" << std::endl;
537 std::cout << "Condition number estimate: " << 1./UMF_Decomposition_Info[UMFPACK_RCOND] << std::endl;
538 std::cout << "Numbers of non-zeroes in decomposition: L: " << UMF_Decomposition_Info[UMFPACK_LNZ] << " U: " << UMF_Decomposition_Info[UMFPACK_UNZ] << std::endl;
539 }
540 if (verbosity_ == 2)
541 {
542 Caller::report_info(UMF_Control,UMF_Decomposition_Info);
543 }
544 }
545
546 void printOnApply(double* UMF_Info)
547 {
548 Caller::report_status(UMF_Control,UMF_Info[UMFPACK_STATUS]);
549 if (verbosity_ > 0)
550 {
551 std::cout << "[UMFPack Solve]" << std::endl;
552 std::cout << "Wallclock Time: " << UMF_Info[UMFPACK_SOLVE_WALLTIME] << " (CPU Time: " << UMF_Info[UMFPACK_SOLVE_TIME] << ")" << std::endl;
553 std::cout << "Flops Taken: " << UMF_Info[UMFPACK_SOLVE_FLOPS] << std::endl;
554 std::cout << "Iterative Refinement steps taken: " << UMF_Info[UMFPACK_IR_TAKEN] << std::endl;
555 std::cout << "Error Estimate: " << UMF_Info[UMFPACK_OMEGA1] << " resp. " << UMF_Info[UMFPACK_OMEGA2] << std::endl;
556 }
557 }
558
559 UMFPackMatrix umfpackMatrix_;
560 bool matrixIsLoaded_;
561 int verbosity_;
562 void *UMF_Symbolic;
563 void *UMF_Numeric;
564 double UMF_Control[UMFPACK_CONTROL];
565 };
566
567 template<typename T, typename A, int n, int m>
568 struct IsDirectSolver<UMFPack<BCRSMatrix<FieldMatrix<T,n,m>,A> > >
569 {
570 enum { value=true};
571 };
572
573 template<typename T, typename A>
574 struct StoresColumnCompressed<UMFPack<BCRSMatrix<T,A> > >
575 {
576 enum { value = true };
577 };
578
579 struct UMFPackCreator {
580 template<class F,class=void> struct isValidBlock : std::false_type{};
581 template<class B> struct isValidBlock<B, std::enable_if_t<std::is_same<typename FieldTraits<B>::real_type,double>::value>> : std::true_type {};
582
583 template<typename TL, typename M>
584 std::shared_ptr<Dune::InverseOperator<typename Dune::TypeListElement<1, TL>::type,
585 typename Dune::TypeListElement<2, TL>::type>>
586 operator() (TL /*tl*/, const M& mat, const Dune::ParameterTree& config,
587 std::enable_if_t<
588 isValidBlock<typename Dune::TypeListElement<1, TL>::type::block_type>::value,int> = 0) const
589 {
590 int verbose = config.get("verbose", 0);
591 return std::make_shared<Dune::UMFPack<M>>(mat,verbose);
592 }
593
594 // second version with SFINAE to validate the template parameters of UMFPack
595 template<typename TL, typename M>
596 std::shared_ptr<Dune::InverseOperator<typename Dune::TypeListElement<1, TL>::type,
597 typename Dune::TypeListElement<2, TL>::type>>
598 operator() (TL /*tl*/, const M& /*mat*/, const Dune::ParameterTree& /*config*/,
599 std::enable_if_t<
600 !isValidBlock<typename Dune::TypeListElement<1, TL>::type::block_type>::value,int> = 0) const
601 {
602 DUNE_THROW(UnsupportedType,
603 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
604 }
605 };
606 DUNE_REGISTER_DIRECT_SOLVER("umfpack",Dune::UMFPackCreator());
607} // end namespace Dune
608
609#endif // HAVE_SUITESPARSE_UMFPACK
610
611#endif //DUNE_ISTL_UMFPACK_HH
Implementation of the BCRSMatrix class.
Inititializer for the ColCompMatrix as needed by OverlappingSchwarz.
Definition: colcompmatrix.hh:255
virtual void setMatrix(const Matrix &mat, const std::set< std::size_t > &mrs)
Initialize data from a given set of matrix rows and columns.
virtual void free()
free allocated space.
size_type N() const
Get the number of rows.
Definition: colcompmatrix.hh:192
size_type M() const
Get the number of columns.
Definition: colcompmatrix.hh:211
Base class for Dune-Exceptions.
Definition: exceptions.hh:94
derive error class from the base class in common
Definition: istlexception.hh:16
Abstract base class for all solvers.
Definition: solver.hh:97
A generic dynamic dense matrix.
Definition: matrix.hh:558
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:214
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_UNUSED_PARAMETER(parm)
A macro to mark intentionally unused function parameters with.
Definition: unused.hh:25
#define DUNE_THROW(E, m)
Definition: exceptions.hh:216
void free()
free allocated space.
Definition: umfpack.hh:491
virtual void apply(domain_type &x, range_type &b, InverseOperatorResult &res)
Apply inverse operator,.
Definition: umfpack.hh:344
virtual SolverCategory::Category category() const
Category of the solver (see SolverCategory::Category)
Definition: umfpack.hh:231
void setMatrix(const Matrix &matrix)
Initialize data from given matrix.
Definition: umfpack.hh:431
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:292
typename Impl::UMFPackVectorChooser< M >::range_type range_type
The type of the range of the solver.
Definition: umfpack.hh:228
UMFPack()
default constructor
Definition: umfpack.hh:274
Dune::ColCompMatrix< Matrix, long int > UMFPackMatrix
The corresponding UMFPack matrix type.
Definition: umfpack.hh:222
void apply(T *x, T *b)
additional apply method with c-arrays in analogy to superlu
Definition: umfpack.hh:386
void setVerbosity(int v)
sets the verbosity level for the UMFPack solver
Definition: umfpack.hh:457
UMFPack(const char *file, int verbose=0)
try loading a decomposition from file
Definition: umfpack.hh:319
void saveDecomposition(const char *file)
saves a decomposition to a file
Definition: umfpack.hh:423
UMFPackMatrix & getInternalMatrix()
Return the column compress matrix from UMFPack.
Definition: umfpack.hh:482
typename Impl::UMFPackVectorChooser< M >::domain_type domain_type
The type of the domain of the solver.
Definition: umfpack.hh:226
UMFPack(const Matrix &matrix, int verbose=0)
Construct a solver object from a matrix.
Definition: umfpack.hh:244
virtual void apply(domain_type &x, range_type &b, double reduction, InverseOperatorResult &res)
apply inverse operator, with given convergence criteria.
Definition: umfpack.hh:375
void setOption(unsigned int option, double value)
Set UMFPack-specific options.
Definition: umfpack.hh:412
M Matrix
The matrix type.
Definition: umfpack.hh:219
ColCompMatrixInitializer< M, long int > MatrixInitializer
Type of an associated initializer class.
Definition: umfpack.hh:224
UMFPack(const Matrix &matrix, int verbose, bool)
Constructor for compatibility with SuperLU standard constructor.
Definition: umfpack.hh:262
void * getFactorization()
Return the matrix factorization.
Definition: umfpack.hh:473
Dune namespace.
Definition: alignedallocator.hh:14
STL namespace.
Implementations of the inverse operator interface.
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
Definition of the DUNE_UNUSED macro for the case that config.h is not available.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 13, 23:29, 2024)