3#ifndef DUNE_ISTL_UMFPACK_HH
4#define DUNE_ISTL_UMFPACK_HH
6#if HAVE_SUITESPARSE_UMFPACK || defined DOXYGEN
16#include<dune/istl/bccsmatrixinitializer.hh>
20#include <dune/istl/solverfactory.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;
46 struct UMFPackMethodChooser
48 static constexpr bool valid = false ;
52 struct UMFPackMethodChooser<double>
54 static constexpr bool valid = true ;
56 template<
typename... A>
57 static void defaults(A... args)
59 umfpack_dl_defaults(args...);
61 template<
typename... A>
62 static void free_numeric(A... args)
64 umfpack_dl_free_numeric(args...);
66 template<
typename... A>
67 static void free_symbolic(A... args)
69 umfpack_dl_free_symbolic(args...);
71 template<
typename... A>
72 static int load_numeric(A... args)
74 return umfpack_dl_load_numeric(args...);
76 template<
typename... A>
77 static void numeric(A... args)
79 umfpack_dl_numeric(args...);
81 template<
typename... A>
82 static void report_info(A... args)
84 umfpack_dl_report_info(args...);
86 template<
typename... A>
87 static void report_status(A... args)
89 umfpack_dl_report_status(args...);
91 template<
typename... A>
92 static int save_numeric(A... args)
94 return umfpack_dl_save_numeric(args...);
96 template<
typename... A>
97 static void solve(A... args)
99 umfpack_dl_solve(args...);
101 template<
typename... A>
102 static void symbolic(A... args)
104 umfpack_dl_symbolic(args...);
109 struct UMFPackMethodChooser<
std::complex<double> >
111 static constexpr bool valid = true ;
113 template<
typename... A>
114 static void defaults(A... args)
116 umfpack_zl_defaults(args...);
118 template<
typename... A>
119 static void free_numeric(A... args)
121 umfpack_zl_free_numeric(args...);
123 template<
typename... A>
124 static void free_symbolic(A... args)
126 umfpack_zl_free_symbolic(args...);
128 template<
typename... A>
129 static int load_numeric(A... args)
131 return umfpack_zl_load_numeric(args...);
133 template<
typename... A>
134 static void numeric(
const long int* cs,
const long int* ri,
const double* val, A... args)
136 umfpack_zl_numeric(cs,ri,val,NULL,args...);
138 template<
typename... A>
139 static void report_info(A... args)
141 umfpack_zl_report_info(args...);
143 template<
typename... A>
144 static void report_status(A... args)
146 umfpack_zl_report_status(args...);
148 template<
typename... A>
149 static int save_numeric(A... args)
151 return umfpack_zl_save_numeric(args...);
153 template<
typename... A>
154 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 const double* cval =
reinterpret_cast<const double*
>(val);
157 umfpack_zl_solve(m,cs,ri,cval,NULL,x,NULL,b,NULL,args...);
159 template<
typename... A>
160 static void symbolic(
long int m,
long int n,
const long int* cs,
const long int* ri,
const double* val, A... args)
162 umfpack_zl_symbolic(m,n,cs,ri,val,NULL,args...);
169 struct UMFPackVectorChooser
172 template<
typename T,
typename A,
int n,
int m>
173 struct UMFPackVectorChooser<BCRSMatrix<FieldMatrix<T,n,m>,A > >
176 using domain_type = BlockVector<
178 typename std::allocator_traits<A>::template rebind_alloc<FieldVector<T,m> > >;
180 using range_type = BlockVector<
182 typename std::allocator_traits<A>::template rebind_alloc<FieldVector<T,n> > >;
185 template<
typename T,
typename A>
186 struct UMFPackVectorChooser<BCRSMatrix<T,A> >
189 using domain_type = BlockVector<T, A>;
191 using range_type = BlockVector<T, A>;
211 typename Impl::UMFPackVectorChooser<M>::domain_type,
212 typename Impl::UMFPackVectorChooser<M>::range_type >
214 using T =
typename M::field_type;
219 using matrix_type = M;
221 typedef ISTL::Impl::BCCSMatrix<typename Matrix::field_type, long int>
UMFPackMatrix;
225 using domain_type =
typename Impl::UMFPackVectorChooser<M>::domain_type;
227 using range_type =
typename Impl::UMFPackVectorChooser<M>::range_type;
232 return SolverCategory::Category::sequential;
246 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
247 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
248 Caller::defaults(UMF_Control);
264 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
265 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
266 Caller::defaults(UMF_Control);
281 :
UMFPack(mat_, config.get<int>(
"verbose", 0))
286 UMFPack() : matrixIsLoaded_(false), verbosity_(0)
289 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
290 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
291 Caller::defaults(UMF_Control);
307 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
308 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
309 Caller::defaults(UMF_Control);
311 int errcode = Caller::load_numeric(&UMF_Numeric,
const_cast<char*
>(file));
312 if ((errcode == UMFPACK_ERROR_out_of_memory) || (errcode == UMFPACK_ERROR_file_IO))
314 matrixIsLoaded_ =
false;
320 matrixIsLoaded_ =
true;
321 std::cout <<
"UMFPack decomposition successfully loaded from " << file << std::endl;
334 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
335 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
336 Caller::defaults(UMF_Control);
337 int errcode = Caller::load_numeric(&UMF_Numeric,
const_cast<char*
>(file));
338 if (errcode == UMFPACK_ERROR_out_of_memory)
340 if (errcode == UMFPACK_ERROR_file_IO)
342 matrixIsLoaded_ =
true;
343 std::cout <<
"UMFPack decomposition successfully loaded from " << file << std::endl;
349 if ((umfpackMatrix_.N() + umfpackMatrix_.M() > 0) || matrixIsLoaded_)
358 if (umfpackMatrix_.N() != b.dim())
360 if (umfpackMatrix_.M() != x.dim())
365 double UMF_Apply_Info[UMFPACK_INFO];
366 Caller::solve(UMFPACK_A,
367 umfpackMatrix_.getColStart(),
368 umfpackMatrix_.getRowIndex(),
369 umfpackMatrix_.getValues(),
370 reinterpret_cast<double*
>(&x[0]),
371 reinterpret_cast<double*
>(&b[0]),
379 res.
elapsed = UMF_Apply_Info[UMFPACK_SOLVE_WALLTIME];
381 printOnApply(UMF_Apply_Info);
399 double UMF_Apply_Info[UMFPACK_INFO];
400 Caller::solve(UMFPACK_A,
401 umfpackMatrix_.getColStart(),
402 umfpackMatrix_.getRowIndex(),
403 umfpackMatrix_.getValues(),
409 printOnApply(UMF_Apply_Info);
425 if (option >= UMFPACK_CONTROL)
428 UMF_Control[option] = value;
436 int errcode = Caller::save_numeric(UMF_Numeric,
const_cast<char*
>(file));
437 if (errcode != UMFPACK_OK)
444 if ((umfpackMatrix_.N() + umfpackMatrix_.M() > 0) || matrixIsLoaded_)
446 if (matrix.N() == 0 or matrix.M() == 0)
449 if (umfpackMatrix_.N() + umfpackMatrix_.M() + umfpackMatrix_.nonzeroes() != 0)
450 umfpackMatrix_.free();
451 umfpackMatrix_.setSize(MatrixDimension<Matrix>::rowdim(matrix),
452 MatrixDimension<Matrix>::coldim(matrix));
453 ISTL::Impl::BCCSMatrixInitializer<Matrix, long int> initializer(umfpackMatrix_);
455 copyToBCCSMatrix(initializer, matrix);
461 void setSubMatrix(
const Matrix& _mat,
const S& rowIndexSet)
463 if ((umfpackMatrix_.N() + umfpackMatrix_.M() > 0) || matrixIsLoaded_)
466 if (umfpackMatrix_.N() + umfpackMatrix_.M() + umfpackMatrix_.nonzeroes() != 0)
467 umfpackMatrix_.free();
469 umfpackMatrix_.setSize(rowIndexSet.size()*MatrixDimension<Matrix>::rowdim(_mat) / _mat.
N(),
470 rowIndexSet.size()*MatrixDimension<Matrix>::coldim(_mat) / _mat.
M());
471 ISTL::Impl::BCCSMatrixInitializer<Matrix, long int> initializer(umfpackMatrix_);
473 copyToBCCSMatrix(initializer, ISTL::Impl::MatrixRowSubset<
Matrix,std::set<std::size_t> >(_mat,rowIndexSet));
490 UMF_Control[UMFPACK_PRL] = 1;
492 UMF_Control[UMFPACK_PRL] = 2;
494 UMF_Control[UMFPACK_PRL] = 4;
512 return umfpackMatrix_;
521 if (!matrixIsLoaded_)
523 Caller::free_symbolic(&UMF_Symbolic);
524 umfpackMatrix_.free();
526 Caller::free_numeric(&UMF_Numeric);
527 matrixIsLoaded_ =
false;
530 const char* name() {
return "UMFPACK"; }
533 typedef typename Dune::UMFPackMethodChooser<T> Caller;
535 template<
class Mat,
class X,
class TM,
class TD,
class T1>
536 friend class SeqOverlappingSchwarz;
537 friend struct SeqOverlappingSchwarzAssemblerHelper<
UMFPack<
Matrix>,true>;
542 double UMF_Decomposition_Info[UMFPACK_INFO];
543 Caller::symbolic(
static_cast<int>(umfpackMatrix_.N()),
544 static_cast<int>(umfpackMatrix_.N()),
545 umfpackMatrix_.getColStart(),
546 umfpackMatrix_.getRowIndex(),
547 reinterpret_cast<double*
>(umfpackMatrix_.getValues()),
550 UMF_Decomposition_Info);
551 Caller::numeric(umfpackMatrix_.getColStart(),
552 umfpackMatrix_.getRowIndex(),
553 reinterpret_cast<double*
>(umfpackMatrix_.getValues()),
557 UMF_Decomposition_Info);
558 Caller::report_status(UMF_Control,UMF_Decomposition_Info[UMFPACK_STATUS]);
561 std::cout <<
"[UMFPack Decomposition]" << std::endl;
562 std::cout <<
"Wallclock Time taken: " << UMF_Decomposition_Info[UMFPACK_NUMERIC_WALLTIME] <<
" (CPU Time: " << UMF_Decomposition_Info[UMFPACK_NUMERIC_TIME] <<
")" << std::endl;
563 std::cout <<
"Flops taken: " << UMF_Decomposition_Info[UMFPACK_FLOPS] << std::endl;
564 std::cout <<
"Peak Memory Usage: " << UMF_Decomposition_Info[UMFPACK_PEAK_MEMORY]*UMF_Decomposition_Info[UMFPACK_SIZE_OF_UNIT] <<
" bytes" << std::endl;
565 std::cout <<
"Condition number estimate: " << 1./UMF_Decomposition_Info[UMFPACK_RCOND] << std::endl;
566 std::cout <<
"Numbers of non-zeroes in decomposition: L: " << UMF_Decomposition_Info[UMFPACK_LNZ] <<
" U: " << UMF_Decomposition_Info[UMFPACK_UNZ] << std::endl;
570 Caller::report_info(UMF_Control,UMF_Decomposition_Info);
574 void printOnApply(
double* UMF_Info)
576 Caller::report_status(UMF_Control,UMF_Info[UMFPACK_STATUS]);
579 std::cout <<
"[UMFPack Solve]" << std::endl;
580 std::cout <<
"Wallclock Time: " << UMF_Info[UMFPACK_SOLVE_WALLTIME] <<
" (CPU Time: " << UMF_Info[UMFPACK_SOLVE_TIME] <<
")" << std::endl;
581 std::cout <<
"Flops Taken: " << UMF_Info[UMFPACK_SOLVE_FLOPS] << std::endl;
582 std::cout <<
"Iterative Refinement steps taken: " << UMF_Info[UMFPACK_IR_TAKEN] << std::endl;
583 std::cout <<
"Error Estimate: " << UMF_Info[UMFPACK_OMEGA1] <<
" resp. " << UMF_Info[UMFPACK_OMEGA2] << std::endl;
588 bool matrixIsLoaded_;
592 double UMF_Control[UMFPACK_CONTROL];
595 template<
typename T,
typename A,
int n,
int m>
596 struct IsDirectSolver<UMFPack<BCRSMatrix<FieldMatrix<T,n,m>,A> > >
601 template<
typename T,
typename A>
602 struct StoresColumnCompressed<UMFPack<BCRSMatrix<T,A> > >
604 enum { value =
true };
607 struct UMFPackCreator {
608 template<
class F,
class=
void>
struct isValidBlock : std::false_type{};
609 template<
class B>
struct isValidBlock<B,
std::enable_if_t<std::is_same<typename FieldTraits<B>::real_type,double>::value>> : std::true_type {};
611 template<
typename TL,
typename M>
612 std::shared_ptr<Dune::InverseOperator<typename Dune::TypeListElement<1, TL>::type,
613 typename Dune::TypeListElement<2, TL>::type>>
616 isValidBlock<
typename Dune::TypeListElement<1, TL>::type::block_type>::value,
int> = 0)
const
618 int verbose = config.
get(
"verbose", 0);
619 return std::make_shared<Dune::UMFPack<M>>(mat,verbose);
623 template<
typename TL,
typename M>
624 std::shared_ptr<Dune::InverseOperator<typename Dune::TypeListElement<1, TL>::type,
625 typename Dune::TypeListElement<2, TL>::type>>
628 !isValidBlock<
typename Dune::TypeListElement<1, TL>::type::block_type>::value,
int> = 0)
const
631 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
634 DUNE_REGISTER_DIRECT_SOLVER(
"umfpack",Dune::UMFPackCreator());
Base class for Dune-Exceptions.
Definition: exceptions.hh:94
derive error class from the base class in common
Definition: istlexception.hh:17
Abstract base class for all solvers.
Definition: solver.hh:97
A generic dynamic dense matrix.
Definition: matrix.hh:559
size_type M() const
Return the number of columns.
Definition: matrix.hh:698
size_type N() const
Return the number of rows.
Definition: matrix.hh:693
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:213
A few common exception classes.
Implementation of the BCRSMatrix class.
Implementations of the inverse operator interface.
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:216
void free()
free allocated space.
Definition: umfpack.hh:519
virtual void apply(domain_type &x, range_type &b, InverseOperatorResult &res)
Apply inverse operator,.
Definition: umfpack.hh:356
virtual SolverCategory::Category category() const
Category of the solver (see SolverCategory::Category)
Definition: umfpack.hh:230
void setMatrix(const Matrix &matrix)
Initialize data from given matrix.
Definition: umfpack.hh:442
UMFPack(const Matrix &mat_, const ParameterTree &config)
Construct a solver object from a matrix.
Definition: umfpack.hh:280
ISTL::Impl::BCCSMatrix< typename Matrix::field_type, long int > UMFPackMatrix
The corresponding UMFPack matrix type.
Definition: umfpack.hh:221
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:304
typename Impl::UMFPackVectorChooser< M >::range_type range_type
The type of the range of the solver.
Definition: umfpack.hh:227
UMFPack()
default constructor
Definition: umfpack.hh:286
ISTL::Impl::BCCSMatrixInitializer< M, long int > MatrixInitializer
Type of an associated initializer class.
Definition: umfpack.hh:223
void apply(T *x, T *b)
additional apply method with c-arrays in analogy to superlu
Definition: umfpack.hh:397
void setVerbosity(int v)
sets the verbosity level for the UMFPack solver
Definition: umfpack.hh:485
UMFPack(const char *file, int verbose=0)
try loading a decomposition from file
Definition: umfpack.hh:331
void saveDecomposition(const char *file)
saves a decomposition to a file
Definition: umfpack.hh:434
UMFPackMatrix & getInternalMatrix()
Return the column compress matrix from UMFPack.
Definition: umfpack.hh:510
typename Impl::UMFPackVectorChooser< M >::domain_type domain_type
The type of the domain of the solver.
Definition: umfpack.hh:225
UMFPack(const Matrix &matrix, int verbose=0)
Construct a solver object from a matrix.
Definition: umfpack.hh:243
virtual void apply(domain_type &x, range_type &b, double reduction, InverseOperatorResult &res)
apply inverse operator, with given convergence criteria.
Definition: umfpack.hh:387
void setOption(unsigned int option, double value)
Set UMFPack-specific options.
Definition: umfpack.hh:423
M Matrix
The matrix type.
Definition: umfpack.hh:218
UMFPack(const Matrix &matrix, int verbose, bool)
Constructor for compatibility with SuperLU standard constructor.
Definition: umfpack.hh:261
void * getFactorization()
Return the matrix factorization.
Definition: umfpack.hh:501
Dune namespace.
Definition: alignedallocator.hh:11
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