Dune Core Modules (2.3.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_UMFPACK_HH
4#define DUNE_UMFPACK_HH
5
6#if HAVE_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
58 template<>
59 struct UMFPackMethodChooser<double>
60 {
61 template<typename... A>
62 static void defaults(A... args)
63 {
64 umfpack_di_defaults(args...);
65 }
66 template<typename... A>
67 static void free_numeric(A... args)
68 {
69 umfpack_di_free_numeric(args...);
70 }
71 template<typename... A>
72 static void free_symbolic(A... args)
73 {
74 umfpack_di_free_symbolic(args...);
75 }
76 template<typename... A>
77 static int load_numeric(A... args)
78 {
79 return umfpack_di_load_numeric(args...);
80 }
81 template<typename... A>
82 static void numeric(A... args)
83 {
84 umfpack_di_numeric(args...);
85 }
86 template<typename... A>
87 static void report_info(A... args)
88 {
89 umfpack_di_report_info(args...);
90 }
91 template<typename... A>
92 static void report_status(A... args)
93 {
94 umfpack_di_report_status(args...);
95 }
96 template<typename... A>
97 static int save_numeric(A... args)
98 {
99 return umfpack_di_save_numeric(args...);
100 }
101 template<typename... A>
102 static void solve(A... args)
103 {
104 umfpack_di_solve(args...);
105 }
106 template<typename... A>
107 static void symbolic(A... args)
108 {
109 umfpack_di_symbolic(args...);
110 }
111 };
112
113 template<>
114 struct UMFPackMethodChooser<std::complex<double> >
115 {
116 template<typename... A>
117 static void defaults(A... args)
118 {
119 umfpack_zi_defaults(args...);
120 }
121 template<typename... A>
122 static void free_numeric(A... args)
123 {
124 umfpack_zi_free_numeric(args...);
125 }
126 template<typename... A>
127 static void free_symbolic(A... args)
128 {
129 umfpack_zi_free_symbolic(args...);
130 }
131 template<typename... A>
132 static int load_numeric(A... args)
133 {
134 return umfpack_zi_load_numeric(args...);
135 }
136 template<typename... A>
137 static void numeric(const int* cs, const int* ri, const double* val, A... args)
138 {
139 umfpack_zi_numeric(cs,ri,val,NULL,args...);
140 }
141 template<typename... A>
142 static void report_info(A... args)
143 {
144 umfpack_zi_report_info(args...);
145 }
146 template<typename... A>
147 static void report_status(A... args)
148 {
149 umfpack_zi_report_status(args...);
150 }
151 template<typename... A>
152 static int save_numeric(A... args)
153 {
154 return umfpack_zi_save_numeric(args...);
155 }
156 template<typename... A>
157 static void solve(int m, const int* cs, const int* ri, std::complex<double>* val, double* x, const double* b,A... args)
158 {
159 const double* cval = reinterpret_cast<const double*>(val);
160 umfpack_zi_solve(m,cs,ri,cval,NULL,x,NULL,b,NULL,args...);
161 }
162 template<typename... A>
163 static void symbolic(int m, int n, const int* cs, const int* ri, const double* val, A... args)
164 {
165 umfpack_zi_symbolic(m,n,cs,ri,val,NULL,args...);
166 }
167 };
168
177 template<typename T, typename A, int n, int m>
178 class UMFPack<BCRSMatrix<FieldMatrix<T,n,m>,A > >
179 : public InverseOperator<
180 BlockVector<FieldVector<T,m>,
181 typename A::template rebind<FieldVector<T,m> >::other>,
182 BlockVector<FieldVector<T,n>,
183 typename A::template rebind<FieldVector<T,n> >::other> >
184 {
185 public:
194 typedef Dune::BlockVector<
196 typename A::template rebind<FieldVector<T,m> >::other> domain_type;
198 typedef Dune::BlockVector<
200 typename A::template rebind<FieldVector<T,n> >::other> range_type;
201
206 UMFPack(const Matrix& mat_, int verbose=0) : mat_is_loaded(false)
207 {
208 //check whether T is a supported type
209 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
210 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
211 Caller::defaults(UMF_Control);
212 setVerbosity(verbose);
213 setMatrix(mat_);
214 }
215
220 UMFPack(const Matrix& mat_, int verbose, bool) : mat_is_loaded(false)
221 {
222 //check whether T is a supported type
223 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
224 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
225 Caller::defaults(UMF_Control);
226 setVerbosity(verbose);
227 setMatrix(mat_);
228 }
229
232 UMFPack() : mat_is_loaded(false), verbose(0)
233 {
234 //check whether T is a supported type
235 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
236 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
237 Caller::defaults(UMF_Control);
238 }
239
249 UMFPack(const Matrix& mat_, const char* file, int verbose=0)
250 {
251 //check whether T is a supported type
252 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
253 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
254 Caller::defaults(UMF_Control);
255 setVerbosity(verbose);
256 int errcode = Caller::load_numeric(&UMF_Numeric, const_cast<char*>(file));
257 if ((errcode == UMFPACK_ERROR_out_of_memory) || (errcode == UMFPACK_ERROR_file_IO))
258 {
259 mat_is_loaded = false;
260 setMatrix(mat_);
261 saveDecomposition(file);
262 }
263 else
264 {
265 mat_is_loaded = true;
266 std::cout << "UMFPack decomposition successfully loaded from " << file << std::endl;
267 }
268 }
269
276 UMFPack(const char* file, int verbose=0)
277 {
278 //check whether T is a supported type
279 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
280 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
281 Caller::defaults(UMF_Control);
282 int errcode = Caller::load_numeric(&UMF_Numeric, const_cast<char*>(file));
283 if (errcode == UMFPACK_ERROR_out_of_memory)
284 DUNE_THROW(Dune::Exception, "ran out of memory while loading UMFPack decomposition");
285 if (errcode == UMFPACK_ERROR_file_IO)
286 DUNE_THROW(Dune::Exception, "IO error while loading UMFPack decomposition");
287 mat_is_loaded = true;
288 std::cout << "UMFPack decomposition successfully loaded from " << file << std::endl;
289 setVerbosity(verbose);
290 }
291
292 virtual ~UMFPack()
293 {
294 if ((mat.N() + mat.M() > 0) || (mat_is_loaded))
295 free();
296 }
297
302 {
303 double UMF_Apply_Info[UMFPACK_INFO];
304 Caller::solve(UMFPACK_A,
305 mat.getColStart(),
306 mat.getRowIndex(),
307 mat.getValues(),
308 reinterpret_cast<double*>(&x[0]),
309 reinterpret_cast<double*>(&b[0]),
310 UMF_Numeric,
311 UMF_Control,
312 UMF_Apply_Info);
313
314 //this is a direct solver
315 res.iterations = 1;
316 res.converged = true;
317
318 printOnApply(UMF_Apply_Info);
319 }
320
324 virtual void apply (domain_type& x, range_type& b, double reduction, InverseOperatorResult& res)
325 {
326 DUNE_UNUSED_PARAMETER(reduction);
327 apply(x,b,res);
328 }
329
335 void apply(T* x, T* b)
336 {
337 double UMF_Apply_Info[UMFPACK_INFO];
338 Caller::solve(UMFPACK_A,
339 mat.getColStart(),
340 mat.getRowIndex(),
341 mat.getValues(),
342 x,
343 b,
344 UMF_Numeric,
345 UMF_Control,
346 UMF_Apply_Info);
347 printOnApply(UMF_Apply_Info);
348 }
349
361 void setOption(unsigned int option, double value)
362 {
363 if (option >= UMFPACK_CONTROL)
364 DUNE_THROW(RangeError, "Requested non-existing UMFPack option");
365
366 UMF_Control[option] = value;
367 }
368
372 void saveDecomposition(const char* file)
373 {
374 int errcode = Caller::save_numeric(UMF_Numeric, const_cast<char*>(file));
375 if (errcode != UMFPACK_OK)
376 DUNE_THROW(Dune::Exception,"IO ERROR while trying to save UMFPack decomposition");
377 }
378
380 void setMatrix(const Matrix& _mat)
381 {
382 if ((mat.N() + mat.M() > 0) || (mat_is_loaded))
383 free();
384 mat = _mat;
385 decompose();
386 }
387
388 template<class S>
389 void setSubMatrix(const Matrix& _mat, const S& rowIndexSet)
390 {
391 if ((mat.N() + mat.M() > 0) || (mat_is_loaded))
392 free();
393 mat.setMatrix(_mat,rowIndexSet);
394 decompose();
395 }
396
404 void setVerbosity(int v)
405 {
406 verbose = v;
407 // set the verbosity level in UMFPack
408 if (verbose == 0)
409 UMF_Control[UMFPACK_PRL] = 1;
410 if (verbose == 1)
411 UMF_Control[UMFPACK_PRL] = 2;
412 if (verbose == 2)
413 UMF_Control[UMFPACK_PRL] = 4;
414 }
415
420 void free()
421 {
422 if (!mat_is_loaded)
423 {
424 Caller::free_symbolic(&UMF_Symbolic);
425 mat.free();
426 }
427 Caller::free_numeric(&UMF_Numeric);
428 mat_is_loaded = false;
429 }
430
431 const char* name() { return "UMFPACK"; }
432
433 private:
434 typedef typename Dune::UMFPackMethodChooser<T> Caller;
435
436 template<class M,class X, class TM, class TD, class T1>
437 friend class SeqOverlappingSchwarz;
438 friend struct SeqOverlappingSchwarzAssemblerHelper<UMFPack<Matrix>,true>;
439
441 void decompose()
442 {
443 double UMF_Decomposition_Info[UMFPACK_INFO];
444 Caller::symbolic(static_cast<int>(mat.N()),
445 static_cast<int>(mat.N()),
446 mat.getColStart(),
447 mat.getRowIndex(),
448 reinterpret_cast<double*>(mat.getValues()),
449 &UMF_Symbolic,
450 UMF_Control,
451 UMF_Decomposition_Info);
452 Caller::numeric(mat.getColStart(),
453 mat.getRowIndex(),
454 reinterpret_cast<double*>(mat.getValues()),
455 UMF_Symbolic,
456 &UMF_Numeric,
457 UMF_Control,
458 UMF_Decomposition_Info);
459 Caller::report_status(UMF_Control,UMF_Decomposition_Info[UMFPACK_STATUS]);
460 if (verbose == 1)
461 {
462 std::cout << "[UMFPack Decomposition]" << std::endl;
463 std::cout << "Wallclock Time taken: " << UMF_Decomposition_Info[UMFPACK_NUMERIC_WALLTIME] << " (CPU Time: " << UMF_Decomposition_Info[UMFPACK_NUMERIC_TIME] << ")" << std::endl;
464 std::cout << "Flops taken: " << UMF_Decomposition_Info[UMFPACK_FLOPS] << std::endl;
465 std::cout << "Peak Memory Usage: " << UMF_Decomposition_Info[UMFPACK_PEAK_MEMORY]*UMF_Decomposition_Info[UMFPACK_SIZE_OF_UNIT] << " bytes" << std::endl;
466 std::cout << "Condition number estimate: " << 1./UMF_Decomposition_Info[UMFPACK_RCOND] << std::endl;
467 std::cout << "Numbers of non-zeroes in decomposition: L: " << UMF_Decomposition_Info[UMFPACK_LNZ] << " U: " << UMF_Decomposition_Info[UMFPACK_UNZ] << std::endl;
468 }
469 if (verbose == 2)
470 {
471 Caller::report_info(UMF_Control,UMF_Decomposition_Info);
472 }
473 }
474
475 void printOnApply(double* UMF_Info)
476 {
477 Caller::report_status(UMF_Control,UMF_Info[UMFPACK_STATUS]);
478 if (verbose > 0)
479 {
480 std::cout << "[UMFPack Solve]" << std::endl;
481 std::cout << "Wallclock Time: " << UMF_Info[UMFPACK_SOLVE_WALLTIME] << " (CPU Time: " << UMF_Info[UMFPACK_SOLVE_TIME] << ")" << std::endl;
482 std::cout << "Flops Taken: " << UMF_Info[UMFPACK_SOLVE_FLOPS] << std::endl;
483 std::cout << "Iterative Refinement steps taken: " << UMF_Info[UMFPACK_IR_TAKEN] << std::endl;
484 std::cout << "Error Estimate: " << UMF_Info[UMFPACK_OMEGA1] << " resp. " << UMF_Info[UMFPACK_OMEGA2] << std::endl;
485 }
486 }
487
488 UMFPackMatrix mat;
489 bool mat_is_loaded;
490 int verbose;
491 void *UMF_Symbolic;
492 void *UMF_Numeric;
493 double UMF_Control[UMFPACK_CONTROL];
494 };
495
496 template<typename T, typename A, int n, int m>
497 struct IsDirectSolver<UMFPack<BCRSMatrix<FieldMatrix<T,n,m>,A> > >
498 {
499 enum { value=true};
500 };
501
502 template<typename T, typename A, int n, int m>
503 struct StoresColumnCompressed<UMFPack<BCRSMatrix<FieldMatrix<T,n,m>,A> > >
504 {
505 enum { value = true };
506 };
507}
508
509#endif //HAVE_UMFPACK
510
511#endif //DUNE_UMFPACK_HH
Implementation of the BCRSMatrix class.
A sparse block matrix with compressed row storage.
Definition: bcrsmatrix.hh:414
A vector of blocks with memory management.
Definition: bvector.hh:254
Base class for Dune-Exceptions.
Definition: exceptions.hh:92
A dense n x m matrix.
Definition: fmatrix.hh:67
vector space out of a tensor product of fields.
Definition: fvector.hh:92
Abstract base class for all solvers.
Definition: solver.hh:79
A generic dynamic dense matrix.
Definition: matrix.hh:25
Default exception class for range errors.
Definition: exceptions.hh:280
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.
#define DUNE_THROW(E, m)
Definition: exceptions.hh:244
void setVerbosity(int v)
sets the verbosity level for the UMFPack solver
Definition: umfpack.hh:404
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:249
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:196
virtual void apply(domain_type &x, range_type &b, InverseOperatorResult &res)
Apply inverse operator,.
Definition: umfpack.hh:301
void setMatrix(const Matrix &_mat)
Initialize data from given matrix.
Definition: umfpack.hh:380
void free()
free allocated space.
Definition: umfpack.hh:420
UMFPack(const Matrix &mat_, int verbose, bool)
Constructor for compatibility with SuperLU standard constructor.
Definition: umfpack.hh:220
ColCompMatrixInitializer< BCRSMatrix< FieldMatrix< T, n, m >, A > > MatrixInitializer
Type of an associated initializer class.
Definition: umfpack.hh:192
void setOption(unsigned int option, double value)
Set UMFPack-specific options.
Definition: umfpack.hh:361
Dune::ColCompMatrix< Matrix > UMFPackMatrix
The corresponding SuperLU Matrix type.
Definition: umfpack.hh:190
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:200
UMFPack(const Matrix &mat_, int verbose=0)
construct a solver object from a BCRSMatrix
Definition: umfpack.hh:206
void apply(T *x, T *b)
additional apply method with c-arrays in analogy to superlu
Definition: umfpack.hh:335
UMFPack()
default constructor
Definition: umfpack.hh:232
UMFPack(const char *file, int verbose=0)
try loading a decomposition from file
Definition: umfpack.hh:276
virtual void apply(domain_type &x, range_type &b, double reduction, InverseOperatorResult &res)
apply inverse operator, with given convergence criteria.
Definition: umfpack.hh:324
void saveDecomposition(const char *file)
saves a decomposition to a file
Definition: umfpack.hh:372
Dune namespace.
Definition: alignment.hh:14
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:32
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 intentional unused function parameters with.
Definition: unused.hh:18
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)