5#ifndef DUNE_ISTL_UMFPACK_HH
6#define DUNE_ISTL_UMFPACK_HH
8#if HAVE_SUITESPARSE_UMFPACK || defined DOXYGEN
18#include<dune/istl/bccsmatrixinitializer.hh>
20#include<dune/istl/foreach.hh>
21#include<dune/istl/multitypeblockmatrix.hh>
22#include<dune/istl/multitypeblockvector.hh>
25#include <dune/istl/solverfactory.hh>
42 template<
class M,
class T,
class TM,
class TD,
class TA>
43 class SeqOverlappingSchwarz;
45 template<
class T,
bool tag>
46 struct SeqOverlappingSchwarzAssemblerHelper;
51 struct UMFPackMethodChooser
53 static constexpr bool valid = false ;
57 struct UMFPackMethodChooser<double>
59 static constexpr bool valid = true ;
61 template<
typename... A>
62 static void defaults(A... args)
64 umfpack_dl_defaults(args...);
66 template<
typename... A>
67 static void free_numeric(A... args)
69 umfpack_dl_free_numeric(args...);
71 template<
typename... A>
72 static void free_symbolic(A... args)
74 umfpack_dl_free_symbolic(args...);
76 template<
typename... A>
77 static int load_numeric(A... args)
79 return umfpack_dl_load_numeric(args...);
81 template<
typename... A>
82 static void numeric(A... args)
84 umfpack_dl_numeric(args...);
86 template<
typename... A>
87 static void report_info(A... args)
89 umfpack_dl_report_info(args...);
91 template<
typename... A>
92 static void report_status(A... args)
94 umfpack_dl_report_status(args...);
96 template<
typename... A>
97 static int save_numeric(A... args)
99 return umfpack_dl_save_numeric(args...);
101 template<
typename... A>
102 static void solve(A... args)
104 umfpack_dl_solve(args...);
106 template<
typename... A>
107 static void symbolic(A... args)
109 umfpack_dl_symbolic(args...);
114 struct UMFPackMethodChooser<
std::complex<double> >
116 static constexpr bool valid = true ;
117 using size_type = SuiteSparse_long;
119 template<
typename... A>
120 static void defaults(A... args)
122 umfpack_zl_defaults(args...);
124 template<
typename... A>
125 static void free_numeric(A... args)
127 umfpack_zl_free_numeric(args...);
129 template<
typename... A>
130 static void free_symbolic(A... args)
132 umfpack_zl_free_symbolic(args...);
134 template<
typename... A>
135 static int load_numeric(A... args)
137 return umfpack_zl_load_numeric(args...);
139 template<
typename... A>
140 static void numeric(
const size_type* cs,
const size_type* ri,
const double* val, A... args)
142 umfpack_zl_numeric(cs,ri,val,NULL,args...);
144 template<
typename... A>
145 static void report_info(A... args)
147 umfpack_zl_report_info(args...);
149 template<
typename... A>
150 static void report_status(A... args)
152 umfpack_zl_report_status(args...);
154 template<
typename... A>
155 static int save_numeric(A... args)
157 return umfpack_zl_save_numeric(args...);
159 template<
typename... A>
160 static void solve(size_type m,
const size_type* cs,
const size_type* ri, std::complex<double>* val,
double* x,
const double* b,A... args)
162 const double* cval =
reinterpret_cast<const double*
>(val);
163 umfpack_zl_solve(m,cs,ri,cval,NULL,x,NULL,b,NULL,args...);
165 template<
typename... A>
166 static void symbolic(size_type m, size_type n,
const size_type* cs,
const size_type* ri,
const double* val, A... args)
168 umfpack_zl_symbolic(m,n,cs,ri,val,NULL,args...);
174 template<
class M,
class =
void>
175 struct UMFPackVectorChooser;
178 template<
class M>
using UMFPackDomainType =
typename UMFPackVectorChooser<M>::domain_type;
181 template<
class M>
using UMFPackRangeType =
typename UMFPackVectorChooser<M>::range_type;
184 struct UMFPackVectorChooser<M,
185 std::enable_if_t<(std::is_same<M,double>::value) || (std::is_same<M,std::complex<double> >::value)>>
187 using domain_type = M;
188 using range_type = M;
191 template<
typename T,
int n,
int m>
192 struct UMFPackVectorChooser<FieldMatrix<T,n,m>,
193 std::enable_if_t<(std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value)>>
196 using domain_type = FieldVector<T,m>;
198 using range_type = FieldVector<T,n>;
201 template<
typename T,
typename A>
202 struct UMFPackVectorChooser<BCRSMatrix<T,A>,
203 std::
void_t<UMFPackDomainType<T>, UMFPackRangeType<T>>>
209 using domain_type = BlockVector<UMFPackDomainType<T>,
typename std::allocator_traits<A>::template rebind_alloc<UMFPackDomainType<T>>>;
211 using range_type = BlockVector<UMFPackRangeType<T>,
typename std::allocator_traits<A>::template rebind_alloc<UMFPackRangeType<T>>>;
215 template<
typename FirstBlock,
typename... Blocks>
216 struct UMFPackVectorChooser<MultiTypeBlockVector<FirstBlock, Blocks...>,
217 std::
void_t<UMFPackDomainType<FirstBlock>, UMFPackRangeType<FirstBlock>, UMFPackDomainType<Blocks>...>>
220 using domain_type = MultiTypeBlockVector<UMFPackDomainType<FirstBlock>, UMFPackDomainType<Blocks>...>;
222 using range_type = UMFPackRangeType<FirstBlock>;
226 template<
typename FirstRow,
typename... Rows>
227 struct UMFPackVectorChooser<MultiTypeBlockMatrix<FirstRow, Rows...>,
228 std::
void_t<UMFPackDomainType<FirstRow>, UMFPackRangeType<FirstRow>, UMFPackRangeType<Rows>...>>
231 using domain_type = UMFPackDomainType<FirstRow>;
233 using range_type = MultiTypeBlockVector< UMFPackRangeType<FirstRow>, UMFPackRangeType<Rows>... >;
259 using T =
typename M::field_type;
262 using size_type = SuiteSparse_long;
266 using matrix_type = M;
268 using UMFPackMatrix = ISTL::Impl::BCCSMatrix<typename Matrix::field_type, size_type>;
279 return SolverCategory::Category::sequential;
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);
311 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
312 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
313 Caller::defaults(UMF_Control);
328 :
UMFPack(mat_, config.
get<int>(
"verbose", 0))
333 UMFPack() : matrixIsLoaded_(false), verbosity_(0)
336 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
337 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
338 Caller::defaults(UMF_Control);
354 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
355 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
356 Caller::defaults(UMF_Control);
358 int errcode = Caller::load_numeric(&UMF_Numeric,
const_cast<char*
>(file));
359 if ((errcode == UMFPACK_ERROR_out_of_memory) || (errcode == UMFPACK_ERROR_file_IO))
361 matrixIsLoaded_ =
false;
367 matrixIsLoaded_ =
true;
368 std::cout <<
"UMFPack decomposition successfully loaded from " << file << std::endl;
381 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
382 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
383 Caller::defaults(UMF_Control);
384 int errcode = Caller::load_numeric(&UMF_Numeric,
const_cast<char*
>(file));
385 if (errcode == UMFPACK_ERROR_out_of_memory)
387 if (errcode == UMFPACK_ERROR_file_IO)
389 matrixIsLoaded_ =
true;
390 std::cout <<
"UMFPack decomposition successfully loaded from " << file << std::endl;
396 if ((umfpackMatrix_.N() + umfpackMatrix_.M() > 0) || matrixIsLoaded_)
405 if (umfpackMatrix_.N() != b.dim())
407 if (umfpackMatrix_.M() != x.dim())
414 std::vector<T> xFlat(x.dim()), bFlat(b.dim());
418 xFlat[ offset ] = entry;
423 bFlat[ offset ] = entry;
426 double UMF_Apply_Info[UMFPACK_INFO];
427 Caller::solve(UMFPACK_A,
428 umfpackMatrix_.getColStart(),
429 umfpackMatrix_.getRowIndex(),
430 umfpackMatrix_.getValues(),
431 reinterpret_cast<double*
>(&xFlat[0]),
432 reinterpret_cast<double*
>(&bFlat[0]),
440 entry = xFlat[offset];
446 res.
elapsed = UMF_Apply_Info[UMFPACK_SOLVE_WALLTIME];
448 printOnApply(UMF_Apply_Info);
468 double UMF_Apply_Info[UMFPACK_INFO];
469 Caller::solve(UMFPACK_A,
470 umfpackMatrix_.getColStart(),
471 umfpackMatrix_.getRowIndex(),
472 umfpackMatrix_.getValues(),
478 printOnApply(UMF_Apply_Info);
494 if (option >= UMFPACK_CONTROL)
497 UMF_Control[option] = value;
505 int errcode = Caller::save_numeric(UMF_Numeric,
const_cast<char*
>(file));
506 if (errcode != UMFPACK_OK)
519 template<
class BitVector = Impl::NoBitVector>
522 if ((umfpackMatrix_.N() + umfpackMatrix_.M() > 0) || matrixIsLoaded_)
524 if (matrix.N() == 0 or matrix.M() == 0)
527 if (umfpackMatrix_.N() + umfpackMatrix_.M() + umfpackMatrix_.nonzeroes() != 0)
528 umfpackMatrix_.free();
530 constexpr bool useBitVector = not std::is_same_v<BitVector,Impl::NoBitVector>;
533 std::vector<bool> flatBitVector;
535 std::vector<size_type> subIndices;
537 [[maybe_unused]]
int numberOfIgnoredDofs = 0;
540 if constexpr ( useBitVector )
543 flatBitVector.resize(flatSize);
547 flatBitVector[ offset ] = entry;
550 numberOfIgnoredDofs++;
556 auto [flatRows,flatCols] =
flatMatrixForEach( matrix, [&](
auto&& ,
auto&& row,
auto&& col){
558 if constexpr ( useBitVector )
559 if ( flatBitVector[row] or flatBitVector[col] )
565 if constexpr ( useBitVector )
570 size_type subIndexCounter = 0;
571 for ( size_type i=0; i<size_type(flatRows); i++ )
572 if ( not flatBitVector[ i ] )
573 subIndices[ i ] = subIndexCounter++;
576 flatRows -= numberOfIgnoredDofs;
577 flatCols -= numberOfIgnoredDofs;
581 umfpackMatrix_.setSize(flatRows,flatCols);
582 umfpackMatrix_.Nnz_ = nonZeros;
585 umfpackMatrix_.colstart =
new size_type[flatCols+1];
586 umfpackMatrix_.rowindex =
new size_type[nonZeros];
587 umfpackMatrix_.values =
new T[nonZeros];
589 for ( size_type i=0; i<size_type(flatCols+1); i++ )
591 umfpackMatrix_.colstart[i] = 0;
599 if constexpr ( useBitVector )
600 if ( flatBitVector[flatRowIndex] or flatBitVector[flatColIndex] )
605 auto colIdx = useBitVector ? subIndices[flatColIndex] : flatColIndex;
607 umfpackMatrix_.colstart[colIdx+1]++;
611 for ( size_type i=0; i<(size_type)flatCols; i++ )
613 umfpackMatrix_.colstart[i+1] += umfpackMatrix_.colstart[i];
617 std::vector<size_type> colPosition(flatCols,0);
620 flatMatrixForEach(matrix, [&](
auto&& entry,
auto&& flatRowIndex,
auto&& flatColIndex)
623 if constexpr ( useBitVector )
624 if ( flatBitVector[flatRowIndex] or flatBitVector[flatColIndex] )
629 auto rowIdx = useBitVector ? subIndices[flatRowIndex] : flatRowIndex;
630 auto colIdx = useBitVector ? subIndices[flatColIndex] : flatColIndex;
633 auto colStart = umfpackMatrix_.colstart[colIdx];
635 auto colPos = colPosition[colIdx];
637 umfpackMatrix_.rowindex[ colStart + colPos ] = rowIdx;
638 umfpackMatrix_.values[ colStart + colPos ] = entry;
640 colPosition[colIdx]++;
651 void setSubMatrix(
const Matrix& _mat,
const S& rowIndexSet)
653 if ((umfpackMatrix_.N() + umfpackMatrix_.M() > 0) || matrixIsLoaded_)
656 if (umfpackMatrix_.N() + umfpackMatrix_.M() + umfpackMatrix_.nonzeroes() != 0)
657 umfpackMatrix_.free();
659 umfpackMatrix_.setSize(rowIndexSet.size()*MatrixDimension<Matrix>::rowdim(_mat) / _mat.N(),
660 rowIndexSet.size()*MatrixDimension<Matrix>::coldim(_mat) / _mat.M());
661 ISTL::Impl::BCCSMatrixInitializer<Matrix, SuiteSparse_long> initializer(umfpackMatrix_);
663 copyToBCCSMatrix(initializer, ISTL::Impl::MatrixRowSubset<
Matrix,std::set<std::size_t> >(_mat,rowIndexSet));
680 UMF_Control[UMFPACK_PRL] = 1;
682 UMF_Control[UMFPACK_PRL] = 2;
684 UMF_Control[UMFPACK_PRL] = 4;
702 return umfpackMatrix_;
711 if (!matrixIsLoaded_)
713 Caller::free_symbolic(&UMF_Symbolic);
714 umfpackMatrix_.free();
716 Caller::free_numeric(&UMF_Numeric);
717 matrixIsLoaded_ =
false;
720 const char* name() {
return "UMFPACK"; }
723 typedef typename Dune::UMFPackMethodChooser<T> Caller;
725 template<
class Mat,
class X,
class TM,
class TD,
class T1>
726 friend class SeqOverlappingSchwarz;
727 friend struct SeqOverlappingSchwarzAssemblerHelper<
UMFPack<
Matrix>,true>;
732 double UMF_Decomposition_Info[UMFPACK_INFO];
733 Caller::symbolic(
static_cast<SuiteSparse_long
>(umfpackMatrix_.N()),
734 static_cast<SuiteSparse_long
>(umfpackMatrix_.N()),
735 umfpackMatrix_.getColStart(),
736 umfpackMatrix_.getRowIndex(),
737 reinterpret_cast<double*
>(umfpackMatrix_.getValues()),
740 UMF_Decomposition_Info);
741 Caller::numeric(umfpackMatrix_.getColStart(),
742 umfpackMatrix_.getRowIndex(),
743 reinterpret_cast<double*
>(umfpackMatrix_.getValues()),
747 UMF_Decomposition_Info);
748 Caller::report_status(UMF_Control,UMF_Decomposition_Info[UMFPACK_STATUS]);
751 std::cout <<
"[UMFPack Decomposition]" << std::endl;
752 std::cout <<
"Wallclock Time taken: " << UMF_Decomposition_Info[UMFPACK_NUMERIC_WALLTIME] <<
" (CPU Time: " << UMF_Decomposition_Info[UMFPACK_NUMERIC_TIME] <<
")" << std::endl;
753 std::cout <<
"Flops taken: " << UMF_Decomposition_Info[UMFPACK_FLOPS] << std::endl;
754 std::cout <<
"Peak Memory Usage: " << UMF_Decomposition_Info[UMFPACK_PEAK_MEMORY]*UMF_Decomposition_Info[UMFPACK_SIZE_OF_UNIT] <<
" bytes" << std::endl;
755 std::cout <<
"Condition number estimate: " << 1./UMF_Decomposition_Info[UMFPACK_RCOND] << std::endl;
756 std::cout <<
"Numbers of non-zeroes in decomposition: L: " << UMF_Decomposition_Info[UMFPACK_LNZ] <<
" U: " << UMF_Decomposition_Info[UMFPACK_UNZ] << std::endl;
760 Caller::report_info(UMF_Control,UMF_Decomposition_Info);
764 void printOnApply(
double* UMF_Info)
766 Caller::report_status(UMF_Control,UMF_Info[UMFPACK_STATUS]);
769 std::cout <<
"[UMFPack Solve]" << std::endl;
770 std::cout <<
"Wallclock Time: " << UMF_Info[UMFPACK_SOLVE_WALLTIME] <<
" (CPU Time: " << UMF_Info[UMFPACK_SOLVE_TIME] <<
")" << std::endl;
771 std::cout <<
"Flops Taken: " << UMF_Info[UMFPACK_SOLVE_FLOPS] << std::endl;
772 std::cout <<
"Iterative Refinement steps taken: " << UMF_Info[UMFPACK_IR_TAKEN] << std::endl;
773 std::cout <<
"Error Estimate: " << UMF_Info[UMFPACK_OMEGA1] <<
" resp. " << UMF_Info[UMFPACK_OMEGA2] << std::endl;
778 bool matrixIsLoaded_;
782 double UMF_Control[UMFPACK_CONTROL];
785 template<
typename T,
typename A,
int n,
int m>
786 struct IsDirectSolver<UMFPack<BCRSMatrix<FieldMatrix<T,n,m>,A> > >
791 template<
typename T,
typename A>
792 struct StoresColumnCompressed<UMFPack<BCRSMatrix<T,A> > >
794 enum { value =
true };
797 namespace UMFPackImpl {
798 template<
class OpTraits,
class=
void>
struct isValidBlock : std::false_type{};
799 template<
class OpTraits>
struct isValidBlock<OpTraits,
801 std::is_same_v<Impl::UMFPackDomainType<typename OpTraits::matrix_type>, typename OpTraits::domain_type>
802 && std::is_same_v<Impl::UMFPackRangeType<typename OpTraits::matrix_type>, typename OpTraits::range_type>
803 && std::is_same_v<typename FieldTraits<typename OpTraits::domain_type::field_type>::real_type, double>
804 && std::is_same_v<typename FieldTraits<typename OpTraits::range_type::field_type>::real_type, double>
805 >> : std::true_type {};
807 DUNE_REGISTER_SOLVER(
"umfpack",
809 -> std::shared_ptr<
typename decltype(opTraits)::solver_type>
811 using OpTraits =
decltype(opTraits);
813 if constexpr (OpTraits::isParallel){
814 if(opTraits.getCommOrThrow(op).communicator().size() > 1)
817 if constexpr (OpTraits::isAssembled){
818 using M =
typename OpTraits::matrix_type;
822 if constexpr (UMFPackImpl::isValidBlock<OpTraits>::value) {
823 const auto& A = opTraits.getAssembledOpOrThrow(op);
824 const M& mat = A->getmat();
825 int verbose = config.
get(
"verbose", 0);
826 return std::make_shared<Dune::UMFPack<M>>(mat,verbose);
830 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
Implementation of the BCRSMatrix class.
Base class for Dune-Exceptions.
Definition: exceptions.hh:96
derive error class from the base class in common
Definition: istlexception.hh:19
Default exception if a function was called while the object is not in a valid state for that function...
Definition: exceptions.hh:281
Abstract base class for all solvers.
Definition: solver.hh:101
Hierarchical structure of string parameters.
Definition: parametertree.hh:37
std::string get(const std::string &key, const std::string &defaultValue) const
get value as string
Definition: parametertree.cc:188
Default exception class for range errors.
Definition: exceptions.hh:254
The UMFPack direct sparse solver.
Definition: umfpack.hh:258
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.
typename Impl::voider< Types... >::type void_t
Is void for all valid input types. The workhorse for C++11 SFINAE-techniques.
Definition: typetraits.hh:40
#define DUNE_THROW(E, m)
Definition: exceptions.hh:218
constexpr auto max
Function object that returns the greater of the given values.
Definition: hybridutilities.hh:484
void free()
free allocated space.
Definition: umfpack.hh:709
UMFPack(const Matrix &mat_, const ParameterTree &config)
Construct a solver object from a matrix.
Definition: umfpack.hh:327
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:351
Impl::UMFPackRangeType< M > range_type
The type of the range of the solver.
Definition: umfpack.hh:274
UMFPack()
default constructor
Definition: umfpack.hh:333
Impl::UMFPackDomainType< M > domain_type
The type of the domain of the solver.
Definition: umfpack.hh:272
void apply(T *x, T *b)
additional apply method with c-arrays in analogy to superlu
Definition: umfpack.hh:466
SolverCategory::Category category() const override
Category of the solver (see SolverCategory::Category)
Definition: umfpack.hh:277
void setVerbosity(int v)
sets the verbosity level for the UMFPack solver
Definition: umfpack.hh:675
UMFPack(const char *file, int verbose=0)
try loading a decomposition from file
Definition: umfpack.hh:378
void apply(domain_type &x, range_type &b, double reduction, InverseOperatorResult &res) override
apply inverse operator, with given convergence criteria.
Definition: umfpack.hh:454
void setMatrix(const Matrix &matrix, const BitVector &bitVector={})
Initialize data from given matrix.
Definition: umfpack.hh:520
void saveDecomposition(const char *file)
saves a decomposition to a file
Definition: umfpack.hh:503
UMFPackMatrix & getInternalMatrix()
Return the column compress matrix from UMFPack.
Definition: umfpack.hh:700
UMFPack(const Matrix &matrix, int verbose=0)
Construct a solver object from a matrix.
Definition: umfpack.hh:290
ISTL::Impl::BCCSMatrixInitializer< M, size_type > MatrixInitializer
Type of an associated initializer class.
Definition: umfpack.hh:270
ISTL::Impl::BCCSMatrix< typename Matrix::field_type, size_type > UMFPackMatrix
The corresponding (scalar) UMFPack matrix type.
Definition: umfpack.hh:268
void setOption(unsigned int option, double value)
Set UMFPack-specific options.
Definition: umfpack.hh:492
M Matrix
The matrix type.
Definition: umfpack.hh:265
void apply(domain_type &x, range_type &b, InverseOperatorResult &res) override
Apply inverse operator,.
Definition: umfpack.hh:403
UMFPack(const Matrix &matrix, int verbose, bool)
Constructor for compatibility with SuperLU standard constructor.
Definition: umfpack.hh:308
void * getFactorization()
Return the matrix factorization.
Definition: umfpack.hh:691
Dune namespace.
Definition: alignedallocator.hh:13
std::pair< std::size_t, std::size_t > flatMatrixForEach(Matrix &&matrix, F &&f, std::size_t rowOffset=0, std::size_t colOffset=0)
Traverse a blocked matrix and call a functor at each scalar entry.
Definition: foreach.hh:132
std::size_t flatVectorForEach(Vector &&vector, F &&f, std::size_t offset=0)
Traverse a blocked vector and call a functor at each scalar entry.
Definition: foreach.hh:95
constexpr auto get(std::integer_sequence< T, II... >, std::integral_constant< std::size_t, pos >={})
Return the entry at position pos of the given sequence.
Definition: integersequence.hh:22
Implementations of the inverse operator interface.
Templates characterizing the type of a solver.
Statistics about the application of an inverse operator.
Definition: solver.hh:50
double elapsed
Elapsed time in seconds.
Definition: solver.hh:84
int iterations
Number of iterations.
Definition: solver.hh:69
bool converged
True if convergence criterion has been met.
Definition: solver.hh:75
Category
Definition: solvercategory.hh:23