3#ifndef DUNE_ISTL_UMFPACK_HH
4#define DUNE_ISTL_UMFPACK_HH
6#if HAVE_SUITESPARSE_UMFPACK || defined DOXYGEN
20#include <dune/istl/solverfactory.hh>
22#include"colcompmatrix.hh"
38 template<
class M,
class T,
class TM,
class TD,
class TA>
39 class SeqOverlappingSchwarz;
41 template<
class T,
bool tag>
42 struct SeqOverlappingSchwarzAssemblerHelper;
47 struct UMFPackMethodChooser
49 static constexpr bool valid = false ;
53 struct UMFPackMethodChooser<double>
55 static constexpr bool valid = true ;
57 template<
typename... A>
58 static void defaults(A... args)
60 umfpack_dl_defaults(args...);
62 template<
typename... A>
63 static void free_numeric(A... args)
65 umfpack_dl_free_numeric(args...);
67 template<
typename... A>
68 static void free_symbolic(A... args)
70 umfpack_dl_free_symbolic(args...);
72 template<
typename... A>
73 static int load_numeric(A... args)
75 return umfpack_dl_load_numeric(args...);
77 template<
typename... A>
78 static void numeric(A... args)
80 umfpack_dl_numeric(args...);
82 template<
typename... A>
83 static void report_info(A... args)
85 umfpack_dl_report_info(args...);
87 template<
typename... A>
88 static void report_status(A... args)
90 umfpack_dl_report_status(args...);
92 template<
typename... A>
93 static int save_numeric(A... args)
95 return umfpack_dl_save_numeric(args...);
97 template<
typename... A>
98 static void solve(A... args)
100 umfpack_dl_solve(args...);
102 template<
typename... A>
103 static void symbolic(A... args)
105 umfpack_dl_symbolic(args...);
110 struct UMFPackMethodChooser<
std::complex<double> >
112 static constexpr bool valid = true ;
114 template<
typename... A>
115 static void defaults(A... args)
117 umfpack_zl_defaults(args...);
119 template<
typename... A>
120 static void free_numeric(A... args)
122 umfpack_zl_free_numeric(args...);
124 template<
typename... A>
125 static void free_symbolic(A... args)
127 umfpack_zl_free_symbolic(args...);
129 template<
typename... A>
130 static int load_numeric(A... args)
132 return umfpack_zl_load_numeric(args...);
134 template<
typename... A>
135 static void numeric(
const long int* cs,
const long int* ri,
const double* val, A... args)
137 umfpack_zl_numeric(cs,ri,val,NULL,args...);
139 template<
typename... A>
140 static void report_info(A... args)
142 umfpack_zl_report_info(args...);
144 template<
typename... A>
145 static void report_status(A... args)
147 umfpack_zl_report_status(args...);
149 template<
typename... A>
150 static int save_numeric(A... args)
152 return umfpack_zl_save_numeric(args...);
154 template<
typename... A>
155 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)
157 const double* cval =
reinterpret_cast<const double*
>(val);
158 umfpack_zl_solve(m,cs,ri,cval,NULL,x,NULL,b,NULL,args...);
160 template<
typename... A>
161 static void symbolic(
long int m,
long int n,
const long int* cs,
const long int* ri,
const double* val, A... args)
163 umfpack_zl_symbolic(m,n,cs,ri,val,NULL,args...);
170 struct UMFPackVectorChooser
173 template<
typename T,
typename A,
int n,
int m>
174 struct UMFPackVectorChooser<BCRSMatrix<FieldMatrix<T,n,m>,A > >
177 using domain_type = BlockVector<
179 typename std::allocator_traits<A>::template rebind_alloc<FieldVector<T,m> > >;
181 using range_type = BlockVector<
183 typename std::allocator_traits<A>::template rebind_alloc<FieldVector<T,n> > >;
186 template<
typename T,
typename A>
187 struct UMFPackVectorChooser<BCRSMatrix<T,A> >
190 using domain_type = BlockVector<T, A>;
192 using range_type = BlockVector<T, A>;
212 typename Impl::UMFPackVectorChooser<M>::domain_type,
213 typename Impl::UMFPackVectorChooser<M>::range_type >
215 using T =
typename M::field_type;
220 using matrix_type = M;
226 using domain_type =
typename Impl::UMFPackVectorChooser<M>::domain_type;
228 using range_type =
typename Impl::UMFPackVectorChooser<M>::range_type;
233 return SolverCategory::Category::sequential;
247 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
248 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
249 Caller::defaults(UMF_Control);
265 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
266 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
267 Caller::defaults(UMF_Control);
274 UMFPack() : matrixIsLoaded_(false), verbosity_(0)
277 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
278 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
279 Caller::defaults(UMF_Control);
295 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
296 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
297 Caller::defaults(UMF_Control);
299 int errcode = Caller::load_numeric(&UMF_Numeric,
const_cast<char*
>(file));
300 if ((errcode == UMFPACK_ERROR_out_of_memory) || (errcode == UMFPACK_ERROR_file_IO))
302 matrixIsLoaded_ =
false;
308 matrixIsLoaded_ =
true;
309 std::cout <<
"UMFPack decomposition successfully loaded from " << file << std::endl;
322 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
323 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
324 Caller::defaults(UMF_Control);
325 int errcode = Caller::load_numeric(&UMF_Numeric,
const_cast<char*
>(file));
326 if (errcode == UMFPACK_ERROR_out_of_memory)
328 if (errcode == UMFPACK_ERROR_file_IO)
330 matrixIsLoaded_ =
true;
331 std::cout <<
"UMFPack decomposition successfully loaded from " << file << std::endl;
337 if ((umfpackMatrix_.
N() + umfpackMatrix_.
M() > 0) || matrixIsLoaded_)
346 if (umfpackMatrix_.
N() != b.dim())
348 if (umfpackMatrix_.
M() != x.dim())
353 double UMF_Apply_Info[UMFPACK_INFO];
354 Caller::solve(UMFPACK_A,
355 umfpackMatrix_.getColStart(),
356 umfpackMatrix_.getRowIndex(),
357 umfpackMatrix_.getValues(),
358 reinterpret_cast<double*
>(&x[0]),
359 reinterpret_cast<double*
>(&b[0]),
367 res.
elapsed = UMF_Apply_Info[UMFPACK_SOLVE_WALLTIME];
369 printOnApply(UMF_Apply_Info);
388 double UMF_Apply_Info[UMFPACK_INFO];
389 Caller::solve(UMFPACK_A,
390 umfpackMatrix_.getColStart(),
391 umfpackMatrix_.getRowIndex(),
392 umfpackMatrix_.getValues(),
398 printOnApply(UMF_Apply_Info);
414 if (option >= UMFPACK_CONTROL)
417 UMF_Control[option] = value;
425 int errcode = Caller::save_numeric(UMF_Numeric,
const_cast<char*
>(file));
426 if (errcode != UMFPACK_OK)
433 if ((umfpackMatrix_.
N() + umfpackMatrix_.
M() > 0) || matrixIsLoaded_)
435 if (matrix.N() == 0 or matrix.M() == 0)
437 umfpackMatrix_ = matrix;
442 void setSubMatrix(
const Matrix& _mat,
const S& rowIndexSet)
444 if ((umfpackMatrix_.
N() + umfpackMatrix_.
M() > 0) || matrixIsLoaded_)
446 umfpackMatrix_.
setMatrix(_mat,rowIndexSet);
462 UMF_Control[UMFPACK_PRL] = 1;
464 UMF_Control[UMFPACK_PRL] = 2;
466 UMF_Control[UMFPACK_PRL] = 4;
484 return umfpackMatrix_;
493 if (!matrixIsLoaded_)
495 Caller::free_symbolic(&UMF_Symbolic);
496 umfpackMatrix_.
free();
498 Caller::free_numeric(&UMF_Numeric);
499 matrixIsLoaded_ =
false;
502 const char* name() {
return "UMFPACK"; }
505 typedef typename Dune::UMFPackMethodChooser<T> Caller;
507 template<
class Mat,
class X,
class TM,
class TD,
class T1>
508 friend class SeqOverlappingSchwarz;
509 friend struct SeqOverlappingSchwarzAssemblerHelper<
UMFPack<
Matrix>,true>;
514 double UMF_Decomposition_Info[UMFPACK_INFO];
515 Caller::symbolic(
static_cast<int>(umfpackMatrix_.
N()),
516 static_cast<int>(umfpackMatrix_.
N()),
517 umfpackMatrix_.getColStart(),
518 umfpackMatrix_.getRowIndex(),
519 reinterpret_cast<double*
>(umfpackMatrix_.getValues()),
522 UMF_Decomposition_Info);
523 Caller::numeric(umfpackMatrix_.getColStart(),
524 umfpackMatrix_.getRowIndex(),
525 reinterpret_cast<double*
>(umfpackMatrix_.getValues()),
529 UMF_Decomposition_Info);
530 Caller::report_status(UMF_Control,UMF_Decomposition_Info[UMFPACK_STATUS]);
533 std::cout <<
"[UMFPack Decomposition]" << std::endl;
534 std::cout <<
"Wallclock Time taken: " << UMF_Decomposition_Info[UMFPACK_NUMERIC_WALLTIME] <<
" (CPU Time: " << UMF_Decomposition_Info[UMFPACK_NUMERIC_TIME] <<
")" << std::endl;
535 std::cout <<
"Flops taken: " << UMF_Decomposition_Info[UMFPACK_FLOPS] << std::endl;
536 std::cout <<
"Peak Memory Usage: " << UMF_Decomposition_Info[UMFPACK_PEAK_MEMORY]*UMF_Decomposition_Info[UMFPACK_SIZE_OF_UNIT] <<
" bytes" << std::endl;
537 std::cout <<
"Condition number estimate: " << 1./UMF_Decomposition_Info[UMFPACK_RCOND] << std::endl;
538 std::cout <<
"Numbers of non-zeroes in decomposition: L: " << UMF_Decomposition_Info[UMFPACK_LNZ] <<
" U: " << UMF_Decomposition_Info[UMFPACK_UNZ] << std::endl;
542 Caller::report_info(UMF_Control,UMF_Decomposition_Info);
546 void printOnApply(
double* UMF_Info)
548 Caller::report_status(UMF_Control,UMF_Info[UMFPACK_STATUS]);
551 std::cout <<
"[UMFPack Solve]" << std::endl;
552 std::cout <<
"Wallclock Time: " << UMF_Info[UMFPACK_SOLVE_WALLTIME] <<
" (CPU Time: " << UMF_Info[UMFPACK_SOLVE_TIME] <<
")" << std::endl;
553 std::cout <<
"Flops Taken: " << UMF_Info[UMFPACK_SOLVE_FLOPS] << std::endl;
554 std::cout <<
"Iterative Refinement steps taken: " << UMF_Info[UMFPACK_IR_TAKEN] << std::endl;
555 std::cout <<
"Error Estimate: " << UMF_Info[UMFPACK_OMEGA1] <<
" resp. " << UMF_Info[UMFPACK_OMEGA2] << std::endl;
560 bool matrixIsLoaded_;
564 double UMF_Control[UMFPACK_CONTROL];
567 template<
typename T,
typename A,
int n,
int m>
568 struct IsDirectSolver<UMFPack<BCRSMatrix<FieldMatrix<T,n,m>,A> > >
573 template<
typename T,
typename A>
574 struct StoresColumnCompressed<UMFPack<BCRSMatrix<T,A> > >
576 enum { value =
true };
579 struct UMFPackCreator {
580 template<
class F,
class=
void>
struct isValidBlock : std::false_type{};
581 template<
class B>
struct isValidBlock<B,
std::enable_if_t<std::is_same<typename FieldTraits<B>::real_type,double>::value>> : std::true_type {};
583 template<
typename TL,
typename M>
584 std::shared_ptr<Dune::InverseOperator<typename Dune::TypeListElement<1, TL>::type,
585 typename Dune::TypeListElement<2, TL>::type>>
588 isValidBlock<
typename Dune::TypeListElement<1, TL>::type::block_type>::value,
int> = 0)
const
590 int verbose = config.
get(
"verbose", 0);
591 return std::make_shared<Dune::UMFPack<M>>(mat,verbose);
595 template<
typename TL,
typename M>
596 std::shared_ptr<Dune::InverseOperator<typename Dune::TypeListElement<1, TL>::type,
597 typename Dune::TypeListElement<2, TL>::type>>
600 !isValidBlock<
typename Dune::TypeListElement<1, TL>::type::block_type>::value,
int> = 0)
const
603 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
606 DUNE_REGISTER_DIRECT_SOLVER(
"umfpack",Dune::UMFPackCreator());
Implementation of the BCRSMatrix class.
Inititializer for the ColCompMatrix as needed by OverlappingSchwarz.
Definition: colcompmatrix.hh:255
virtual void setMatrix(const Matrix &mat, const std::set< std::size_t > &mrs)
Initialize data from a given set of matrix rows and columns.
virtual void free()
free allocated space.
size_type N() const
Get the number of rows.
Definition: colcompmatrix.hh:192
size_type M() const
Get the number of columns.
Definition: colcompmatrix.hh:211
Base class for Dune-Exceptions.
Definition: exceptions.hh:94
derive error class from the base class in common
Definition: istlexception.hh:16
Abstract base class for all solvers.
Definition: solver.hh:97
A generic dynamic dense matrix.
Definition: matrix.hh:558
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:214
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_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
void free()
free allocated space.
Definition: umfpack.hh:491
virtual void apply(domain_type &x, range_type &b, InverseOperatorResult &res)
Apply inverse operator,.
Definition: umfpack.hh:344
virtual SolverCategory::Category category() const
Category of the solver (see SolverCategory::Category)
Definition: umfpack.hh:231
void setMatrix(const Matrix &matrix)
Initialize data from given matrix.
Definition: umfpack.hh:431
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:292
typename Impl::UMFPackVectorChooser< M >::range_type range_type
The type of the range of the solver.
Definition: umfpack.hh:228
UMFPack()
default constructor
Definition: umfpack.hh:274
Dune::ColCompMatrix< Matrix, long int > UMFPackMatrix
The corresponding UMFPack matrix type.
Definition: umfpack.hh:222
void apply(T *x, T *b)
additional apply method with c-arrays in analogy to superlu
Definition: umfpack.hh:386
void setVerbosity(int v)
sets the verbosity level for the UMFPack solver
Definition: umfpack.hh:457
UMFPack(const char *file, int verbose=0)
try loading a decomposition from file
Definition: umfpack.hh:319
void saveDecomposition(const char *file)
saves a decomposition to a file
Definition: umfpack.hh:423
UMFPackMatrix & getInternalMatrix()
Return the column compress matrix from UMFPack.
Definition: umfpack.hh:482
typename Impl::UMFPackVectorChooser< M >::domain_type domain_type
The type of the domain of the solver.
Definition: umfpack.hh:226
UMFPack(const Matrix &matrix, int verbose=0)
Construct a solver object from a matrix.
Definition: umfpack.hh:244
virtual void apply(domain_type &x, range_type &b, double reduction, InverseOperatorResult &res)
apply inverse operator, with given convergence criteria.
Definition: umfpack.hh:375
void setOption(unsigned int option, double value)
Set UMFPack-specific options.
Definition: umfpack.hh:412
M Matrix
The matrix type.
Definition: umfpack.hh:219
ColCompMatrixInitializer< M, long int > MatrixInitializer
Type of an associated initializer class.
Definition: umfpack.hh:224
UMFPack(const Matrix &matrix, int verbose, bool)
Constructor for compatibility with SuperLU standard constructor.
Definition: umfpack.hh:262
void * getFactorization()
Return the matrix factorization.
Definition: umfpack.hh:473
Dune namespace.
Definition: alignedallocator.hh:14
Implementations of the inverse operator interface.
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
Definition of the DUNE_UNUSED macro for the case that config.h is not available.