Dune Core Modules (2.6.0)

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
21#include"colcompmatrix.hh"
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
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
215 {
216 return SolverCategory::Category::sequential;
217 }
218
227 UMFPack(const Matrix& matrix, int verbose=0) : matrixIsLoaded_(false)
228 {
229 //check whether T is a supported type
230 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
231 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
232 Caller::defaults(UMF_Control);
233 setVerbosity(verbose);
234 setMatrix(matrix);
235 }
236
245 UMFPack(const Matrix& matrix, int verbose, bool) : 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
257 UMFPack() : matrixIsLoaded_(false), verbosity_(0)
258 {
259 //check whether T is a supported type
260 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
261 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
262 Caller::defaults(UMF_Control);
263 }
264
275 UMFPack(const Matrix& mat_, const char* file, int verbose=0)
276 {
277 //check whether T is a supported type
278 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
279 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
280 Caller::defaults(UMF_Control);
281 setVerbosity(verbose);
282 int errcode = Caller::load_numeric(&UMF_Numeric, const_cast<char*>(file));
283 if ((errcode == UMFPACK_ERROR_out_of_memory) || (errcode == UMFPACK_ERROR_file_IO))
284 {
285 matrixIsLoaded_ = false;
286 setMatrix(mat_);
287 saveDecomposition(file);
288 }
289 else
290 {
291 matrixIsLoaded_ = true;
292 std::cout << "UMFPack decomposition successfully loaded from " << file << std::endl;
293 }
294 }
295
302 UMFPack(const char* file, int verbose=0)
303 {
304 //check whether T is a supported type
305 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
306 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
307 Caller::defaults(UMF_Control);
308 int errcode = Caller::load_numeric(&UMF_Numeric, const_cast<char*>(file));
309 if (errcode == UMFPACK_ERROR_out_of_memory)
310 DUNE_THROW(Dune::Exception, "ran out of memory while loading UMFPack decomposition");
311 if (errcode == UMFPACK_ERROR_file_IO)
312 DUNE_THROW(Dune::Exception, "IO error while loading UMFPack decomposition");
313 matrixIsLoaded_ = true;
314 std::cout << "UMFPack decomposition successfully loaded from " << file << std::endl;
315 setVerbosity(verbose);
316 }
317
318 virtual ~UMFPack()
319 {
320 if ((umfpackMatrix_.N() + umfpackMatrix_.M() > 0) || matrixIsLoaded_)
321 free();
322 }
323
328 {
329 if (umfpackMatrix_.N() != b.dim())
330 DUNE_THROW(Dune::ISTLError, "Size of right-hand-side vector b does not match the number of matrix rows!");
331 if (umfpackMatrix_.M() != x.dim())
332 DUNE_THROW(Dune::ISTLError, "Size of solution vector x does not match the number of matrix columns!");
333
334 double UMF_Apply_Info[UMFPACK_INFO];
335 Caller::solve(UMFPACK_A,
336 umfpackMatrix_.getColStart(),
337 umfpackMatrix_.getRowIndex(),
338 umfpackMatrix_.getValues(),
339 reinterpret_cast<double*>(&x[0]),
340 reinterpret_cast<double*>(&b[0]),
341 UMF_Numeric,
342 UMF_Control,
343 UMF_Apply_Info);
344
345 //this is a direct solver
346 res.iterations = 1;
347 res.converged = true;
348 res.elapsed = UMF_Apply_Info[UMFPACK_SOLVE_WALLTIME];
349
350 printOnApply(UMF_Apply_Info);
351 }
352
356 virtual void apply (domain_type& x, range_type& b, double reduction, InverseOperatorResult& res)
357 {
358 DUNE_UNUSED_PARAMETER(reduction);
359 apply(x,b,res);
360 }
361
367 void apply(T* x, T* b)
368 {
369 double UMF_Apply_Info[UMFPACK_INFO];
370 Caller::solve(UMFPACK_A,
371 umfpackMatrix_.getColStart(),
372 umfpackMatrix_.getRowIndex(),
373 umfpackMatrix_.getValues(),
374 x,
375 b,
376 UMF_Numeric,
377 UMF_Control,
378 UMF_Apply_Info);
379 printOnApply(UMF_Apply_Info);
380 }
381
393 void setOption(unsigned int option, double value)
394 {
395 if (option >= UMFPACK_CONTROL)
396 DUNE_THROW(RangeError, "Requested non-existing UMFPack option");
397
398 UMF_Control[option] = value;
399 }
400
404 void saveDecomposition(const char* file)
405 {
406 int errcode = Caller::save_numeric(UMF_Numeric, const_cast<char*>(file));
407 if (errcode != UMFPACK_OK)
408 DUNE_THROW(Dune::Exception,"IO ERROR while trying to save UMFPack decomposition");
409 }
410
412 void setMatrix(const Matrix& matrix)
413 {
414 if ((umfpackMatrix_.N() + umfpackMatrix_.M() > 0) || matrixIsLoaded_)
415 free();
416 umfpackMatrix_ = matrix;
417 decompose();
418 }
419
420 template<class S>
421 void setSubMatrix(const Matrix& _mat, const S& rowIndexSet)
422 {
423 if ((umfpackMatrix_.N() + umfpackMatrix_.M() > 0) || matrixIsLoaded_)
424 free();
425 umfpackMatrix_.setMatrix(_mat,rowIndexSet);
426 decompose();
427 }
428
436 void setVerbosity(int v)
437 {
438 verbosity_ = v;
439 // set the verbosity level in UMFPack
440 if (verbosity_ == 0)
441 UMF_Control[UMFPACK_PRL] = 1;
442 if (verbosity_ == 1)
443 UMF_Control[UMFPACK_PRL] = 2;
444 if (verbosity_ == 2)
445 UMF_Control[UMFPACK_PRL] = 4;
446 }
447
453 {
454 return UMF_Numeric;
455 }
456
462 {
463 return umfpackMatrix_;
464 }
465
470 void free()
471 {
472 if (!matrixIsLoaded_)
473 {
474 Caller::free_symbolic(&UMF_Symbolic);
475 umfpackMatrix_.free();
476 }
477 Caller::free_numeric(&UMF_Numeric);
478 matrixIsLoaded_ = false;
479 }
480
481 const char* name() { return "UMFPACK"; }
482
483 private:
484 typedef typename Dune::UMFPackMethodChooser<T> Caller;
485
486 template<class M,class X, class TM, class TD, class T1>
487 friend class SeqOverlappingSchwarz;
488 friend struct SeqOverlappingSchwarzAssemblerHelper<UMFPack<Matrix>,true>;
489
491 void decompose()
492 {
493 double UMF_Decomposition_Info[UMFPACK_INFO];
494 Caller::symbolic(static_cast<int>(umfpackMatrix_.N()),
495 static_cast<int>(umfpackMatrix_.N()),
496 umfpackMatrix_.getColStart(),
497 umfpackMatrix_.getRowIndex(),
498 reinterpret_cast<double*>(umfpackMatrix_.getValues()),
499 &UMF_Symbolic,
500 UMF_Control,
501 UMF_Decomposition_Info);
502 Caller::numeric(umfpackMatrix_.getColStart(),
503 umfpackMatrix_.getRowIndex(),
504 reinterpret_cast<double*>(umfpackMatrix_.getValues()),
505 UMF_Symbolic,
506 &UMF_Numeric,
507 UMF_Control,
508 UMF_Decomposition_Info);
509 Caller::report_status(UMF_Control,UMF_Decomposition_Info[UMFPACK_STATUS]);
510 if (verbosity_ == 1)
511 {
512 std::cout << "[UMFPack Decomposition]" << std::endl;
513 std::cout << "Wallclock Time taken: " << UMF_Decomposition_Info[UMFPACK_NUMERIC_WALLTIME] << " (CPU Time: " << UMF_Decomposition_Info[UMFPACK_NUMERIC_TIME] << ")" << std::endl;
514 std::cout << "Flops taken: " << UMF_Decomposition_Info[UMFPACK_FLOPS] << std::endl;
515 std::cout << "Peak Memory Usage: " << UMF_Decomposition_Info[UMFPACK_PEAK_MEMORY]*UMF_Decomposition_Info[UMFPACK_SIZE_OF_UNIT] << " bytes" << std::endl;
516 std::cout << "Condition number estimate: " << 1./UMF_Decomposition_Info[UMFPACK_RCOND] << std::endl;
517 std::cout << "Numbers of non-zeroes in decomposition: L: " << UMF_Decomposition_Info[UMFPACK_LNZ] << " U: " << UMF_Decomposition_Info[UMFPACK_UNZ] << std::endl;
518 }
519 if (verbosity_ == 2)
520 {
521 Caller::report_info(UMF_Control,UMF_Decomposition_Info);
522 }
523 }
524
525 void printOnApply(double* UMF_Info)
526 {
527 Caller::report_status(UMF_Control,UMF_Info[UMFPACK_STATUS]);
528 if (verbosity_ > 0)
529 {
530 std::cout << "[UMFPack Solve]" << std::endl;
531 std::cout << "Wallclock Time: " << UMF_Info[UMFPACK_SOLVE_WALLTIME] << " (CPU Time: " << UMF_Info[UMFPACK_SOLVE_TIME] << ")" << std::endl;
532 std::cout << "Flops Taken: " << UMF_Info[UMFPACK_SOLVE_FLOPS] << std::endl;
533 std::cout << "Iterative Refinement steps taken: " << UMF_Info[UMFPACK_IR_TAKEN] << std::endl;
534 std::cout << "Error Estimate: " << UMF_Info[UMFPACK_OMEGA1] << " resp. " << UMF_Info[UMFPACK_OMEGA2] << std::endl;
535 }
536 }
537
538 UMFPackMatrix umfpackMatrix_;
539 bool matrixIsLoaded_;
540 int verbosity_;
541 void *UMF_Symbolic;
542 void *UMF_Numeric;
543 double UMF_Control[UMFPACK_CONTROL];
544 };
545
546 template<typename T, typename A, int n, int m>
547 struct IsDirectSolver<UMFPack<BCRSMatrix<FieldMatrix<T,n,m>,A> > >
548 {
549 enum { value=true};
550 };
551
552 template<typename T, typename A, int n, int m>
553 struct StoresColumnCompressed<UMFPack<BCRSMatrix<FieldMatrix<T,n,m>,A> > >
554 {
555 enum { value = true };
556 };
557}
558
559#endif // HAVE_SUITESPARSE_UMFPACK
560
561#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:317
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:91
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
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.
decltype(auto) apply(F &&f, ArgTuple &&args)
Apply function with arguments given as tuple.
Definition: apply.hh:58
#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
UMFPack(const Matrix &matrix, int verbose, bool)
Constructor for compatibility with SuperLU standard constructor.
Definition: umfpack.hh:245
void setVerbosity(int v)
sets the verbosity level for the UMFPack solver
Definition: umfpack.hh:436
void * getFactorization()
Return the matrix factorization.
Definition: umfpack.hh:452
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:275
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
void setMatrix(const Matrix &matrix)
Initialize data from given matrix.
Definition: umfpack.hh:412
virtual void apply(domain_type &x, range_type &b, InverseOperatorResult &res)
Apply inverse operator,.
Definition: umfpack.hh:327
void free()
free allocated space.
Definition: umfpack.hh:470
ColCompMatrixInitializer< BCRSMatrix< FieldMatrix< T, n, m >, A > > MatrixInitializer
Type of an associated initializer class.
Definition: umfpack.hh:203
UMFPackMatrix & getInternalMatrix()
Return the column compress matrix from UMFPack.
Definition: umfpack.hh:461
UMFPack(const Matrix &matrix, int verbose=0)
Construct a solver object from a BCRSMatrix.
Definition: umfpack.hh:227
void setOption(unsigned int option, double value)
Set UMFPack-specific options.
Definition: umfpack.hh:393
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:367
UMFPack()
default constructor
Definition: umfpack.hh:257
UMFPack(const char *file, int verbose=0)
try loading a decomposition from file
Definition: umfpack.hh:302
virtual SolverCategory::Category category() const
Category of the solver (see SolverCategory::Category)
Definition: umfpack.hh:214
virtual void apply(domain_type &x, range_type &b, double reduction, InverseOperatorResult &res)
apply inverse operator, with given convergence criteria.
Definition: umfpack.hh:356
void saveDecomposition(const char *file)
saves a decomposition to a file
Definition: umfpack.hh:404
Dune namespace.
Definition: alignedallocator.hh:10
STL namespace.
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:41
double elapsed
Elapsed time in seconds.
Definition: solver.hh:74
int iterations
Number of iterations.
Definition: solver.hh:59
bool converged
True if convergence criterion has been met.
Definition: solver.hh:65
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 (Dec 26, 23:30, 2024)