3 #ifndef DUNE_ISTL_UMFPACK_HH
4 #define DUNE_ISTL_UMFPACK_HH
8 #if ENABLE_UMFPACK && not defined(ENABLE_SUITESPARSE)
9 #define ENABLE_SUITESPARSE ENABLE_UMFPACK
12 #if HAVE_UMFPACK || defined DOXYGEN
19 #include<dune/common/exceptions.hh>
20 #include<dune/common/fmatrix.hh>
21 #include<dune/common/fvector.hh>
22 #include<dune/common/unused.hh>
43 template<
class M,
class T,
class TM,
class TD,
class TA>
44 class SeqOverlappingSchwarz;
46 template<
class T,
bool tag>
47 struct SeqOverlappingSchwarzAssemblerHelper;
54 template<
class Matrix>
67 template<
typename...
A>
70 umfpack_di_defaults(args...);
72 template<
typename...
A>
75 umfpack_di_free_numeric(args...);
77 template<
typename...
A>
80 umfpack_di_free_symbolic(args...);
82 template<
typename...
A>
85 return umfpack_di_load_numeric(args...);
87 template<
typename...
A>
90 umfpack_di_numeric(args...);
92 template<
typename...
A>
95 umfpack_di_report_info(args...);
97 template<
typename...
A>
100 umfpack_di_report_status(args...);
102 template<
typename...
A>
105 return umfpack_di_save_numeric(args...);
107 template<
typename...
A>
110 umfpack_di_solve(args...);
112 template<
typename...
A>
115 umfpack_di_symbolic(args...);
122 template<
typename...
A>
125 umfpack_zi_defaults(args...);
127 template<
typename...
A>
130 umfpack_zi_free_numeric(args...);
132 template<
typename...
A>
135 umfpack_zi_free_symbolic(args...);
137 template<
typename...
A>
140 return umfpack_zi_load_numeric(args...);
142 template<
typename...
A>
143 static void numeric(
const int* cs,
const int* ri,
const double* val,
A... args)
145 umfpack_zi_numeric(cs,ri,val,NULL,args...);
147 template<
typename...
A>
150 umfpack_zi_report_info(args...);
152 template<
typename...
A>
155 umfpack_zi_report_status(args...);
157 template<
typename...
A>
160 return umfpack_zi_save_numeric(args...);
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)
165 const double* cval =
reinterpret_cast<const double*
>(val);
166 umfpack_zi_solve(m,cs,ri,cval,NULL,x,NULL,b,NULL,args...);
168 template<
typename...
A>
169 static void symbolic(
int m,
int n,
const int* cs,
const int* ri,
const double* val,
A... args)
171 umfpack_zi_symbolic(m,n,cs,ri,val,NULL,args...);
188 template<
typename T,
typename A,
int n,
int m>
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> >
207 typename A::template rebind<FieldVector<T,m> >::other>
domain_type;
211 typename A::template rebind<FieldVector<T,n> >::other>
range_type;
224 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
225 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
226 Caller::defaults(UMF_Control);
227 setVerbosity(verbose);
242 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
243 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
244 Caller::defaults(UMF_Control);
245 setVerbosity(verbose);
251 UMFPack() : matrixIsLoaded_(false), verbosity_(0)
254 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
255 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
256 Caller::defaults(UMF_Control);
272 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
273 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
274 Caller::defaults(UMF_Control);
275 setVerbosity(verbose);
276 int errcode = Caller::load_numeric(&UMF_Numeric, const_cast<char*>(file));
277 if ((errcode == UMFPACK_ERROR_out_of_memory) || (errcode == UMFPACK_ERROR_file_IO))
279 matrixIsLoaded_ =
false;
281 saveDecomposition(file);
285 matrixIsLoaded_ =
true;
286 std::cout <<
"UMFPack decomposition successfully loaded from " << file << std::endl;
299 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
300 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
301 Caller::defaults(UMF_Control);
302 int errcode = Caller::load_numeric(&UMF_Numeric, const_cast<char*>(file));
303 if (errcode == UMFPACK_ERROR_out_of_memory)
304 DUNE_THROW(Dune::Exception,
"ran out of memory while loading UMFPack decomposition");
305 if (errcode == UMFPACK_ERROR_file_IO)
306 DUNE_THROW(Dune::Exception,
"IO error while loading UMFPack decomposition");
307 matrixIsLoaded_ =
true;
308 std::cout <<
"UMFPack decomposition successfully loaded from " << file << std::endl;
309 setVerbosity(verbose);
314 if ((umfpackMatrix_.N() + umfpackMatrix_.M() > 0) || (matrixIsLoaded_))
323 if (umfpackMatrix_.N() != b.
dim())
324 DUNE_THROW(
Dune::ISTLError,
"Size of right-hand-side vector b does not match the number of matrix rows!");
325 if (umfpackMatrix_.M() != x.
dim())
326 DUNE_THROW(
Dune::ISTLError,
"Size of solution vector x does not match the number of matrix columns!");
328 double UMF_Apply_Info[UMFPACK_INFO];
329 Caller::solve(UMFPACK_A,
330 umfpackMatrix_.getColStart(),
331 umfpackMatrix_.getRowIndex(),
332 umfpackMatrix_.getValues(),
333 reinterpret_cast<double*
>(&x[0]),
334 reinterpret_cast<double*>(&b[0]),
342 res.
elapsed = UMF_Apply_Info[UMFPACK_SOLVE_WALLTIME];
344 printOnApply(UMF_Apply_Info);
352 DUNE_UNUSED_PARAMETER(reduction);
363 double UMF_Apply_Info[UMFPACK_INFO];
364 Caller::solve(UMFPACK_A,
365 umfpackMatrix_.getColStart(),
366 umfpackMatrix_.getRowIndex(),
367 umfpackMatrix_.getValues(),
373 printOnApply(UMF_Apply_Info);
389 if (option >= UMFPACK_CONTROL)
390 DUNE_THROW(RangeError,
"Requested non-existing UMFPack option");
392 UMF_Control[option] = value;
400 int errcode = Caller::save_numeric(UMF_Numeric, const_cast<char*>(file));
401 if (errcode != UMFPACK_OK)
402 DUNE_THROW(Dune::Exception,
"IO ERROR while trying to save UMFPack decomposition");
408 if ((umfpackMatrix_.N() + umfpackMatrix_.M() > 0) || (matrixIsLoaded_))
410 umfpackMatrix_ = matrix;
417 if ((umfpackMatrix_.N() + umfpackMatrix_.M() > 0) || (matrixIsLoaded_))
419 umfpackMatrix_.setMatrix(_mat,rowIndexSet);
435 UMF_Control[UMFPACK_PRL] = 1;
437 UMF_Control[UMFPACK_PRL] = 2;
439 UMF_Control[UMFPACK_PRL] = 4;
448 if (!matrixIsLoaded_)
450 Caller::free_symbolic(&UMF_Symbolic);
451 umfpackMatrix_.free();
453 Caller::free_numeric(&UMF_Numeric);
454 matrixIsLoaded_ =
false;
457 const char*
name() {
return "UMFPACK"; }
462 template<
class M,
class X,
class TM,
class TD,
class T1>
469 double UMF_Decomposition_Info[UMFPACK_INFO];
470 Caller::symbolic(static_cast<int>(umfpackMatrix_.N()),
471 static_cast<int>(umfpackMatrix_.N()),
472 umfpackMatrix_.getColStart(),
473 umfpackMatrix_.getRowIndex(),
474 reinterpret_cast<double*
>(umfpackMatrix_.getValues()),
477 UMF_Decomposition_Info);
478 Caller::numeric(umfpackMatrix_.getColStart(),
479 umfpackMatrix_.getRowIndex(),
480 reinterpret_cast<double*
>(umfpackMatrix_.getValues()),
484 UMF_Decomposition_Info);
485 Caller::report_status(UMF_Control,UMF_Decomposition_Info[UMFPACK_STATUS]);
488 std::cout <<
"[UMFPack Decomposition]" << std::endl;
489 std::cout <<
"Wallclock Time taken: " << UMF_Decomposition_Info[UMFPACK_NUMERIC_WALLTIME] <<
" (CPU Time: " << UMF_Decomposition_Info[UMFPACK_NUMERIC_TIME] <<
")" << std::endl;
490 std::cout <<
"Flops taken: " << UMF_Decomposition_Info[UMFPACK_FLOPS] << std::endl;
491 std::cout <<
"Peak Memory Usage: " << UMF_Decomposition_Info[UMFPACK_PEAK_MEMORY]*UMF_Decomposition_Info[UMFPACK_SIZE_OF_UNIT] <<
" bytes" << std::endl;
492 std::cout <<
"Condition number estimate: " << 1./UMF_Decomposition_Info[UMFPACK_RCOND] << std::endl;
493 std::cout <<
"Numbers of non-zeroes in decomposition: L: " << UMF_Decomposition_Info[UMFPACK_LNZ] <<
" U: " << UMF_Decomposition_Info[UMFPACK_UNZ] << std::endl;
497 Caller::report_info(UMF_Control,UMF_Decomposition_Info);
501 void printOnApply(
double* UMF_Info)
503 Caller::report_status(UMF_Control,UMF_Info[UMFPACK_STATUS]);
506 std::cout <<
"[UMFPack Solve]" << std::endl;
507 std::cout <<
"Wallclock Time: " << UMF_Info[UMFPACK_SOLVE_WALLTIME] <<
" (CPU Time: " << UMF_Info[UMFPACK_SOLVE_TIME] <<
")" << std::endl;
508 std::cout <<
"Flops Taken: " << UMF_Info[UMFPACK_SOLVE_FLOPS] << std::endl;
509 std::cout <<
"Iterative Refinement steps taken: " << UMF_Info[UMFPACK_IR_TAKEN] << std::endl;
510 std::cout <<
"Error Estimate: " << UMF_Info[UMFPACK_OMEGA1] <<
" resp. " << UMF_Info[UMFPACK_OMEGA2] << std::endl;
514 UMFPackMatrix& getInternalMatrix() {
return umfpackMatrix_; }
516 UMFPackMatrix umfpackMatrix_;
517 bool matrixIsLoaded_;
521 double UMF_Control[UMFPACK_CONTROL];
524 template<
typename T,
typename A,
int n,
int m>
530 template<
typename T,
typename A,
int n,
int m>
537 #endif //HAVE_UMFPACK
539 #endif //DUNE_ISTL_UMFPACK_HH
Dune::BCRSMatrix< FieldMatrix< T, n, m >, A > Matrix
The matrix type.
Definition: umfpack.hh:198
Dune::BCRSMatrix< FieldMatrix< T, n, m >, A > matrix_type
Definition: umfpack.hh:199
static void free_numeric(A...args)
Definition: umfpack.hh:128
size_type dim() const
dimension of the vector space
Definition: bvector.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:269
ColCompMatrixInitializer< BCRSMatrix< FieldMatrix< T, n, m >, A > > MatrixInitializer
Type of an associated initializer class.
Definition: umfpack.hh:203
Definition: umfpack.hh:61
static void solve(int m, const int *cs, const int *ri, std::complex< double > *val, double *x, const double *b, A...args)
Definition: umfpack.hh:163
UMFPack(const Matrix &matrix, int verbose=0)
Construct a solver object from a BCRSMatrix.
Definition: umfpack.hh:221
derive error class from the base class in common
Definition: istlexception.hh:16
static void free_numeric(A...args)
Definition: umfpack.hh:73
Inititializer for the ColCompMatrix as needed by OverlappingSchwarz.
Definition: colcompmatrix.hh:153
A vector of blocks with memory management.
Definition: bvector.hh:252
double elapsed
Elapsed time in seconds.
Definition: solver.hh:62
int iterations
Number of iterations.
Definition: solver.hh:50
Sequential overlapping Schwarz preconditioner.
Definition: colcompmatrix.hh:157
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
static void report_status(A...args)
Definition: umfpack.hh:98
Matrix & A
Definition: matrixmatrix.hh:216
A sparse block matrix with compressed row storage.
Definition: bcrsmatrix.hh:412
static int save_numeric(A...args)
Definition: umfpack.hh:158
static void defaults(A...args)
Definition: umfpack.hh:123
Dune::ColCompMatrix< Matrix > UMFPackMatrix
The corresponding SuperLU Matrix type.
Definition: umfpack.hh:201
static void free_symbolic(A...args)
Definition: umfpack.hh:133
Templates characterizing the type of a solver.
Definition: solvertype.hh:27
void apply(T *x, T *b)
additional apply method with c-arrays in analogy to superlu
Definition: umfpack.hh:361
const char * name()
Definition: umfpack.hh:457
static int load_numeric(A...args)
Definition: umfpack.hh:138
Whether this is a direct solver.
Definition: solvertype.hh:22
void saveDecomposition(const char *file)
saves a decomposition to a file
Definition: umfpack.hh:398
void free()
free allocated space.
Definition: umfpack.hh:446
void setVerbosity(int v)
sets the verbosity level for the UMFPack solver
Definition: umfpack.hh:430
UMFPack(const Matrix &matrix, int verbose, bool)
Constructor for compatibility with SuperLU standard constructor.
Definition: umfpack.hh:239
Definition: overlappingschwarz.hh:690
static void numeric(A...args)
Definition: umfpack.hh:88
virtual ~UMFPack()
Definition: umfpack.hh:312
void setMatrix(const Matrix &matrix)
Initialize data from given matrix.
Definition: umfpack.hh:406
void setOption(unsigned int option, double value)
Set UMFPack-specific options.
Definition: umfpack.hh:387
Definition: matrixutils.hh:25
virtual void apply(domain_type &x, range_type &b, InverseOperatorResult &res)
Apply inverse operator,.
Definition: umfpack.hh:321
virtual void apply(domain_type &x, range_type &b, double reduction, InverseOperatorResult &res)
apply inverse operator, with given convergence criteria.
Definition: umfpack.hh:350
Abstract base class for all solvers.
Definition: solver.hh:79
whether the solver internally uses column compressed storage
Definition: solvertype.hh:34
void setSubMatrix(const Matrix &_mat, const S &rowIndexSet)
Definition: umfpack.hh:415
Definition: solvertype.hh:13
static void free_symbolic(A...args)
Definition: umfpack.hh:78
static void symbolic(A...args)
Definition: umfpack.hh:113
Use the UMFPack package to directly solve linear systems – empty default class.
Definition: umfpack.hh:55
static void solve(A...args)
Definition: umfpack.hh:108
Implementation of the BCRSMatrix class.
static int load_numeric(A...args)
Definition: umfpack.hh:83
static void numeric(const int *cs, const int *ri, const double *val, A...args)
Definition: umfpack.hh:143
static void symbolic(int m, int n, const int *cs, const int *ri, const double *val, A...args)
Definition: umfpack.hh:169
UMFPack()
default constructor
Definition: umfpack.hh:251
Implementations of the inverse operator interface.
static void report_info(A...args)
Definition: umfpack.hh:93
static int save_numeric(A...args)
Definition: umfpack.hh:103
static void defaults(A...args)
Definition: umfpack.hh:68
static void report_status(A...args)
Definition: umfpack.hh:153
bool converged
True if convergence criterion has been met.
Definition: solver.hh:56
UMFPack(const char *file, int verbose=0)
try loading a decomposition from file
Definition: umfpack.hh:296
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
Statistics about the application of an inverse operator.
Definition: solver.hh:31
static void report_info(A...args)
Definition: umfpack.hh:148