3 #ifndef DUNE_ISTL_UMFPACK_HH
4 #define DUNE_ISTL_UMFPACK_HH
6 #if HAVE_UMFPACK || defined DOXYGEN
13 #include<dune/common/exceptions.hh>
14 #include<dune/common/fmatrix.hh>
15 #include<dune/common/fvector.hh>
16 #include<dune/common/unused.hh>
37 template<
class M,
class T,
class TM,
class TD,
class TA>
38 class SeqOverlappingSchwarz;
40 template<
class T,
bool tag>
41 struct SeqOverlappingSchwarzAssemblerHelper;
48 template<
class Matrix>
61 template<
typename...
A>
64 umfpack_di_defaults(args...);
66 template<
typename...
A>
69 umfpack_di_free_numeric(args...);
71 template<
typename...
A>
74 umfpack_di_free_symbolic(args...);
76 template<
typename...
A>
79 return umfpack_di_load_numeric(args...);
81 template<
typename...
A>
84 umfpack_di_numeric(args...);
86 template<
typename...
A>
89 umfpack_di_report_info(args...);
91 template<
typename...
A>
94 umfpack_di_report_status(args...);
96 template<
typename...
A>
99 return umfpack_di_save_numeric(args...);
101 template<
typename...
A>
104 umfpack_di_solve(args...);
106 template<
typename...
A>
109 umfpack_di_symbolic(args...);
116 template<
typename...
A>
119 umfpack_zi_defaults(args...);
121 template<
typename...
A>
124 umfpack_zi_free_numeric(args...);
126 template<
typename...
A>
129 umfpack_zi_free_symbolic(args...);
131 template<
typename...
A>
134 return umfpack_zi_load_numeric(args...);
136 template<
typename...
A>
137 static void numeric(
const int* cs,
const int* ri,
const double* val,
A... args)
139 umfpack_zi_numeric(cs,ri,val,NULL,args...);
141 template<
typename...
A>
144 umfpack_zi_report_info(args...);
146 template<
typename...
A>
149 umfpack_zi_report_status(args...);
151 template<
typename...
A>
154 return umfpack_zi_save_numeric(args...);
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)
159 const double* cval =
reinterpret_cast<const double*
>(val);
160 umfpack_zi_solve(m,cs,ri,cval,NULL,x,NULL,b,NULL,args...);
162 template<
typename...
A>
163 static void symbolic(
int m,
int n,
const int* cs,
const int* ri,
const double* val,
A... args)
165 umfpack_zi_symbolic(m,n,cs,ri,val,NULL,args...);
182 template<
typename T,
typename A,
int n,
int m>
185 BlockVector<FieldVector<T,m>,
186 typename A::template rebind<FieldVector<T,m> >::other>,
187 BlockVector<FieldVector<T,n>,
188 typename A::template rebind<FieldVector<T,n> >::other> >
201 typename A::template rebind<FieldVector<T,m> >::other>
domain_type;
205 typename A::template rebind<FieldVector<T,n> >::other>
range_type;
215 UMFPack(
const Matrix& matrix,
int verbose=0) : matrixIsLoaded_(false)
218 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
219 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
220 Caller::defaults(UMF_Control);
221 setVerbosity(verbose);
233 UMFPack(
const Matrix& matrix,
int verbose,
bool) : matrixIsLoaded_(false)
236 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
237 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
238 Caller::defaults(UMF_Control);
239 setVerbosity(verbose);
245 UMFPack() : matrixIsLoaded_(false), verbose(0)
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);
263 UMFPack(
const Matrix& mat_,
const char* file,
int verbose=0)
266 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
267 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
268 Caller::defaults(UMF_Control);
269 setVerbosity(verbose);
270 int errcode = Caller::load_numeric(&UMF_Numeric, const_cast<char*>(file));
271 if ((errcode == UMFPACK_ERROR_out_of_memory) || (errcode == UMFPACK_ERROR_file_IO))
273 matrixIsLoaded_ =
false;
275 saveDecomposition(file);
279 matrixIsLoaded_ =
true;
280 std::cout <<
"UMFPack decomposition successfully loaded from " << file << std::endl;
293 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
294 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
295 Caller::defaults(UMF_Control);
296 int errcode = Caller::load_numeric(&UMF_Numeric, const_cast<char*>(file));
297 if (errcode == UMFPACK_ERROR_out_of_memory)
298 DUNE_THROW(Dune::Exception,
"ran out of memory while loading UMFPack decomposition");
299 if (errcode == UMFPACK_ERROR_file_IO)
300 DUNE_THROW(Dune::Exception,
"IO error while loading UMFPack decomposition");
301 matrixIsLoaded_ =
true;
302 std::cout <<
"UMFPack decomposition successfully loaded from " << file << std::endl;
303 setVerbosity(verbose);
308 if ((umfpackMatrix_.N() + umfpackMatrix_.M() > 0) || (matrixIsLoaded_))
317 double UMF_Apply_Info[UMFPACK_INFO];
318 Caller::solve(UMFPACK_A,
319 umfpackMatrix_.getColStart(),
320 umfpackMatrix_.getRowIndex(),
321 umfpackMatrix_.getValues(),
322 reinterpret_cast<double*
>(&x[0]),
323 reinterpret_cast<double*>(&b[0]),
331 res.
elapsed = UMF_Apply_Info[UMFPACK_SOLVE_WALLTIME];
333 printOnApply(UMF_Apply_Info);
341 DUNE_UNUSED_PARAMETER(reduction);
352 double UMF_Apply_Info[UMFPACK_INFO];
353 Caller::solve(UMFPACK_A,
354 umfpackMatrix_.getColStart(),
355 umfpackMatrix_.getRowIndex(),
356 umfpackMatrix_.getValues(),
362 printOnApply(UMF_Apply_Info);
378 if (option >= UMFPACK_CONTROL)
379 DUNE_THROW(RangeError,
"Requested non-existing UMFPack option");
381 UMF_Control[option] = value;
389 int errcode = Caller::save_numeric(UMF_Numeric, const_cast<char*>(file));
390 if (errcode != UMFPACK_OK)
391 DUNE_THROW(Dune::Exception,
"IO ERROR while trying to save UMFPack decomposition");
397 if ((umfpackMatrix_.N() + umfpackMatrix_.M() > 0) || (matrixIsLoaded_))
399 umfpackMatrix_ = matrix;
406 if ((umfpackMatrix_.N() + umfpackMatrix_.M() > 0) || (matrixIsLoaded_))
408 umfpackMatrix_.setMatrix(_mat,rowIndexSet);
424 UMF_Control[UMFPACK_PRL] = 1;
426 UMF_Control[UMFPACK_PRL] = 2;
428 UMF_Control[UMFPACK_PRL] = 4;
437 if (!matrixIsLoaded_)
439 Caller::free_symbolic(&UMF_Symbolic);
440 umfpackMatrix_.free();
442 Caller::free_numeric(&UMF_Numeric);
443 matrixIsLoaded_ =
false;
446 const char*
name() {
return "UMFPACK"; }
451 template<
class M,
class X,
class TM,
class TD,
class T1>
458 double UMF_Decomposition_Info[UMFPACK_INFO];
459 Caller::symbolic(static_cast<int>(umfpackMatrix_.N()),
460 static_cast<int>(umfpackMatrix_.N()),
461 umfpackMatrix_.getColStart(),
462 umfpackMatrix_.getRowIndex(),
463 reinterpret_cast<double*
>(umfpackMatrix_.getValues()),
466 UMF_Decomposition_Info);
467 Caller::numeric(umfpackMatrix_.getColStart(),
468 umfpackMatrix_.getRowIndex(),
469 reinterpret_cast<double*
>(umfpackMatrix_.getValues()),
473 UMF_Decomposition_Info);
474 Caller::report_status(UMF_Control,UMF_Decomposition_Info[UMFPACK_STATUS]);
477 std::cout <<
"[UMFPack Decomposition]" << std::endl;
478 std::cout <<
"Wallclock Time taken: " << UMF_Decomposition_Info[UMFPACK_NUMERIC_WALLTIME] <<
" (CPU Time: " << UMF_Decomposition_Info[UMFPACK_NUMERIC_TIME] <<
")" << std::endl;
479 std::cout <<
"Flops taken: " << UMF_Decomposition_Info[UMFPACK_FLOPS] << std::endl;
480 std::cout <<
"Peak Memory Usage: " << UMF_Decomposition_Info[UMFPACK_PEAK_MEMORY]*UMF_Decomposition_Info[UMFPACK_SIZE_OF_UNIT] <<
" bytes" << std::endl;
481 std::cout <<
"Condition number estimate: " << 1./UMF_Decomposition_Info[UMFPACK_RCOND] << std::endl;
482 std::cout <<
"Numbers of non-zeroes in decomposition: L: " << UMF_Decomposition_Info[UMFPACK_LNZ] <<
" U: " << UMF_Decomposition_Info[UMFPACK_UNZ] << std::endl;
486 Caller::report_info(UMF_Control,UMF_Decomposition_Info);
490 void printOnApply(
double* UMF_Info)
492 Caller::report_status(UMF_Control,UMF_Info[UMFPACK_STATUS]);
495 std::cout <<
"[UMFPack Solve]" << std::endl;
496 std::cout <<
"Wallclock Time: " << UMF_Info[UMFPACK_SOLVE_WALLTIME] <<
" (CPU Time: " << UMF_Info[UMFPACK_SOLVE_TIME] <<
")" << std::endl;
497 std::cout <<
"Flops Taken: " << UMF_Info[UMFPACK_SOLVE_FLOPS] << std::endl;
498 std::cout <<
"Iterative Refinement steps taken: " << UMF_Info[UMFPACK_IR_TAKEN] << std::endl;
499 std::cout <<
"Error Estimate: " << UMF_Info[UMFPACK_OMEGA1] <<
" resp. " << UMF_Info[UMFPACK_OMEGA2] << std::endl;
503 UMFPackMatrix& getInternalMatrix() {
return umfpackMatrix_; }
505 UMFPackMatrix umfpackMatrix_;
506 bool matrixIsLoaded_;
510 double UMF_Control[UMFPACK_CONTROL];
513 template<
typename T,
typename A,
int n,
int m>
519 template<
typename T,
typename A,
int n,
int m>
526 #endif //HAVE_UMFPACK
528 #endif //DUNE_ISTL_UMFPACK_HH
whether the solver internally uses column compressed storage
Definition: solvertype.hh:34
Use the UMFPack package to directly solve linear systems – empty default class.
Definition: umfpack.hh:49
Dune::BCRSMatrix< FieldMatrix< T, n, m >, A > matrix_type
Definition: umfpack.hh:193
Implementation of the BCRSMatrix class.
virtual void apply(domain_type &x, range_type &b, double reduction, InverseOperatorResult &res)
apply inverse operator, with given convergence criteria.
Definition: umfpack.hh:339
void free()
free allocated space.
Definition: umfpack.hh:435
Definition: basearray.hh:19
static void free_numeric(A...args)
Definition: umfpack.hh:122
static void report_info(A...args)
Definition: umfpack.hh:142
ColCompMatrixInitializer< BCRSMatrix< FieldMatrix< T, n, m >, A > > MatrixInitializer
Type of an associated initializer class.
Definition: umfpack.hh:197
static void free_symbolic(A...args)
Definition: umfpack.hh:127
UMFPack(const Matrix &matrix, int verbose, bool)
Constructor for compatibility with SuperLU standard constructor.
Definition: umfpack.hh:233
Inititializer for the ColCompMatrix as needed by OverlappingSchwarz.
Definition: colcompmatrix.hh:153
UMFPack(const char *file, int verbose=0)
try loading a decomposition from file
Definition: umfpack.hh:290
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
static void report_status(A...args)
Definition: umfpack.hh:92
static void defaults(A...args)
Definition: umfpack.hh:62
A sparse block matrix with compressed row storage.
Definition: bcrsmatrix.hh:412
Definition: overlappingschwarz.hh:690
Definition: umfpack.hh:55
static int save_numeric(A...args)
Definition: umfpack.hh:152
Matrix & A
Definition: matrixmatrix.hh:216
UMFPack(const Matrix &matrix, int verbose=0)
Construct a solver object from a BCRSMatrix.
Definition: umfpack.hh:215
static void numeric(const int *cs, const int *ri, const double *val, A...args)
Definition: umfpack.hh:137
static void report_status(A...args)
Definition: umfpack.hh:147
static void solve(A...args)
Definition: umfpack.hh:102
Whether this is a direct solver.
Definition: solvertype.hh:22
void saveDecomposition(const char *file)
saves a decomposition to a file
Definition: umfpack.hh:387
static void free_symbolic(A...args)
Definition: umfpack.hh:72
void setSubMatrix(const Matrix &_mat, const S &rowIndexSet)
Definition: umfpack.hh:404
virtual void apply(domain_type &x, range_type &b, InverseOperatorResult &res)
Apply inverse operator,.
Definition: umfpack.hh:315
void apply(T *x, T *b)
additional apply method with c-arrays in analogy to superlu
Definition: umfpack.hh:350
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:201
Definition: matrixutils.hh:24
Sequential overlapping Schwarz preconditioner.
Definition: colcompmatrix.hh:157
static int load_numeric(A...args)
Definition: umfpack.hh:77
void setOption(unsigned int option, double value)
Set UMFPack-specific options.
Definition: umfpack.hh:376
Abstract base class for all solvers.
Definition: solver.hh:79
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:157
void setVerbosity(int v)
sets the verbosity level for the UMFPack solver
Definition: umfpack.hh:419
Definition: solvertype.hh:13
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:205
static void report_info(A...args)
Definition: umfpack.hh:87
Dune::BCRSMatrix< FieldMatrix< T, n, m >, A > Matrix
The matrix type.
Definition: umfpack.hh:192
static void symbolic(int m, int n, const int *cs, const int *ri, const double *val, A...args)
Definition: umfpack.hh:163
Dune::ColCompMatrix< Matrix > UMFPackMatrix
The corresponding SuperLU Matrix type.
Definition: umfpack.hh:195
Templates characterizing the type of a solver.
Definition: solvertype.hh:27
void setMatrix(const Matrix &matrix)
Initialize data from given matrix.
Definition: umfpack.hh:395
virtual ~UMFPack()
Definition: umfpack.hh:306
UMFPack()
default constructor
Definition: umfpack.hh:245
static int load_numeric(A...args)
Definition: umfpack.hh:132
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:263
static void numeric(A...args)
Definition: umfpack.hh:82
bool converged
True if convergence criterion has been met.
Definition: solver.hh:56
static void symbolic(A...args)
Definition: umfpack.hh:107
Statistics about the application of an inverse operator.
Definition: solver.hh:31
static void defaults(A...args)
Definition: umfpack.hh:117
static void free_numeric(A...args)
Definition: umfpack.hh:67
static int save_numeric(A...args)
Definition: umfpack.hh:97
const char * name()
Definition: umfpack.hh:446
Implementations of the inverse operator interface.