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...);
175 struct UMFPackVectorChooser
178 template<
typename T,
typename A,
int n,
int m>
179 struct UMFPackVectorChooser<BCRSMatrix<FieldMatrix<T,n,m>,A > >
182 using domain_type = BlockVector<
184 typename std::allocator_traits<A>::template rebind_alloc<FieldVector<T,m> > >;
186 using range_type = BlockVector<
188 typename std::allocator_traits<A>::template rebind_alloc<FieldVector<T,n> > >;
191 template<
typename T,
typename A>
192 struct UMFPackVectorChooser<BCRSMatrix<T,A> >
195 using domain_type = BlockVector<T, A>;
197 using range_type = BlockVector<T, A>;
201 template<
typename FirstBlock,
typename... Blocks>
202 struct UMFPackVectorChooser<MultiTypeBlockVector<FirstBlock, Blocks...> >
205 using domain_type = MultiTypeBlockVector<typename UMFPackVectorChooser<FirstBlock>::domain_type,
typename UMFPackVectorChooser<Blocks>::domain_type... >;
207 using range_type =
typename UMFPackVectorChooser<FirstBlock>::range_type;
211 template<
typename FirstRow,
typename... Rows>
212 struct UMFPackVectorChooser<MultiTypeBlockMatrix<FirstRow, Rows...> >
215 using domain_type =
typename UMFPackVectorChooser<FirstRow>::domain_type;
217 using range_type = MultiTypeBlockVector< typename UMFPackVectorChooser<FirstRow>::range_type,
typename UMFPackVectorChooser<Rows>::range_type... >;
243 typename Impl::UMFPackVectorChooser<M>::domain_type,
244 typename Impl::UMFPackVectorChooser<M>::range_type >
246 using T =
typename M::field_type;
250 using size_type = SuiteSparse_long;
254 using matrix_type = M;
256 using UMFPackMatrix = ISTL::Impl::BCCSMatrix<typename Matrix::field_type, size_type>;
260 using domain_type =
typename Impl::UMFPackVectorChooser<M>::domain_type;
262 using range_type =
typename Impl::UMFPackVectorChooser<M>::range_type;
267 return SolverCategory::Category::sequential;
281 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
282 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
283 Caller::defaults(UMF_Control);
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);
316 :
UMFPack(mat_, config.
get<int>(
"verbose", 0))
321 UMFPack() : matrixIsLoaded_(false), verbosity_(0)
324 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
325 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
326 Caller::defaults(UMF_Control);
342 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
343 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
344 Caller::defaults(UMF_Control);
346 int errcode = Caller::load_numeric(&UMF_Numeric,
const_cast<char*
>(file));
347 if ((errcode == UMFPACK_ERROR_out_of_memory) || (errcode == UMFPACK_ERROR_file_IO))
349 matrixIsLoaded_ =
false;
355 matrixIsLoaded_ =
true;
356 std::cout <<
"UMFPack decomposition successfully loaded from " << file << std::endl;
369 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
370 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
371 Caller::defaults(UMF_Control);
372 int errcode = Caller::load_numeric(&UMF_Numeric,
const_cast<char*
>(file));
373 if (errcode == UMFPACK_ERROR_out_of_memory)
375 if (errcode == UMFPACK_ERROR_file_IO)
377 matrixIsLoaded_ =
true;
378 std::cout <<
"UMFPack decomposition successfully loaded from " << file << std::endl;
384 if ((umfpackMatrix_.N() + umfpackMatrix_.M() > 0) || matrixIsLoaded_)
393 if (umfpackMatrix_.N() != b.dim())
395 if (umfpackMatrix_.M() != x.dim())
402 std::vector<T> xFlat(x.dim()), bFlat(b.dim());
406 xFlat[ offset ] = entry;
411 bFlat[ offset ] = entry;
414 double UMF_Apply_Info[UMFPACK_INFO];
415 Caller::solve(UMFPACK_A,
416 umfpackMatrix_.getColStart(),
417 umfpackMatrix_.getRowIndex(),
418 umfpackMatrix_.getValues(),
419 reinterpret_cast<double*
>(&xFlat[0]),
420 reinterpret_cast<double*
>(&bFlat[0]),
428 entry = xFlat[offset];
434 res.
elapsed = UMF_Apply_Info[UMFPACK_SOLVE_WALLTIME];
436 printOnApply(UMF_Apply_Info);
456 double UMF_Apply_Info[UMFPACK_INFO];
457 Caller::solve(UMFPACK_A,
458 umfpackMatrix_.getColStart(),
459 umfpackMatrix_.getRowIndex(),
460 umfpackMatrix_.getValues(),
466 printOnApply(UMF_Apply_Info);
482 if (option >= UMFPACK_CONTROL)
485 UMF_Control[option] = value;
493 int errcode = Caller::save_numeric(UMF_Numeric,
const_cast<char*
>(file));
494 if (errcode != UMFPACK_OK)
507 template<
class BitVector = Impl::NoBitVector>
510 if ((umfpackMatrix_.N() + umfpackMatrix_.M() > 0) || matrixIsLoaded_)
512 if (matrix.N() == 0 or matrix.M() == 0)
515 if (umfpackMatrix_.N() + umfpackMatrix_.M() + umfpackMatrix_.nonzeroes() != 0)
516 umfpackMatrix_.free();
518 constexpr
bool useBitVector = not std::is_same_v<BitVector,Impl::NoBitVector>;
521 std::vector<bool> flatBitVector;
523 std::vector<size_type> subIndices;
525 int numberOfIgnoredDofs = 0;
528 if constexpr ( useBitVector )
531 flatBitVector.resize(flatSize);
535 flatBitVector[ offset ] = entry;
538 numberOfIgnoredDofs++;
544 auto [flatRows,flatCols] =
flatMatrixForEach( matrix, [&](
auto&& ,
auto&& row,
auto&& col){
546 if constexpr ( useBitVector )
547 if ( flatBitVector[row] or flatBitVector[col] )
553 if constexpr ( useBitVector )
558 size_type subIndexCounter = 0;
559 for ( size_type i=0; i<size_type(flatRows); i++ )
560 if ( not flatBitVector[ i ] )
561 subIndices[ i ] = subIndexCounter++;
564 flatRows -= numberOfIgnoredDofs;
565 flatCols -= numberOfIgnoredDofs;
569 umfpackMatrix_.setSize(flatRows,flatCols);
570 umfpackMatrix_.Nnz_ = nonZeros;
573 umfpackMatrix_.colstart =
new size_type[flatCols+1];
574 umfpackMatrix_.rowindex =
new size_type[nonZeros];
575 umfpackMatrix_.values =
new T[nonZeros];
577 for ( size_type i=0; i<size_type(flatCols+1); i++ )
579 umfpackMatrix_.colstart[i] = 0;
587 if constexpr ( useBitVector )
588 if ( flatBitVector[flatRowIndex] or flatBitVector[flatColIndex] )
593 auto colIdx = useBitVector ? subIndices[flatColIndex] : flatColIndex;
595 umfpackMatrix_.colstart[colIdx+1]++;
599 for ( size_type i=0; i<(size_type)flatCols; i++ )
601 umfpackMatrix_.colstart[i+1] += umfpackMatrix_.colstart[i];
605 std::vector<size_type> colPosition(flatCols,0);
608 flatMatrixForEach(matrix, [&](
auto&& entry,
auto&& flatRowIndex,
auto&& flatColIndex)
611 if constexpr ( useBitVector )
612 if ( flatBitVector[flatRowIndex] or flatBitVector[flatColIndex] )
617 auto rowIdx = useBitVector ? subIndices[flatRowIndex] : flatRowIndex;
618 auto colIdx = useBitVector ? subIndices[flatColIndex] : flatColIndex;
621 auto colStart = umfpackMatrix_.colstart[colIdx];
623 auto colPos = colPosition[colIdx];
625 umfpackMatrix_.rowindex[ colStart + colPos ] = rowIdx;
626 umfpackMatrix_.values[ colStart + colPos ] = entry;
628 colPosition[colIdx]++;
639 void setSubMatrix(
const Matrix& _mat,
const S& rowIndexSet)
641 if ((umfpackMatrix_.N() + umfpackMatrix_.M() > 0) || matrixIsLoaded_)
644 if (umfpackMatrix_.N() + umfpackMatrix_.M() + umfpackMatrix_.nonzeroes() != 0)
645 umfpackMatrix_.free();
647 umfpackMatrix_.setSize(rowIndexSet.size()*MatrixDimension<Matrix>::rowdim(_mat) / _mat.N(),
648 rowIndexSet.size()*MatrixDimension<Matrix>::coldim(_mat) / _mat.M());
649 ISTL::Impl::BCCSMatrixInitializer<Matrix, SuiteSparse_long> initializer(umfpackMatrix_);
651 copyToBCCSMatrix(initializer, ISTL::Impl::MatrixRowSubset<
Matrix,std::set<std::size_t> >(_mat,rowIndexSet));
668 UMF_Control[UMFPACK_PRL] = 1;
670 UMF_Control[UMFPACK_PRL] = 2;
672 UMF_Control[UMFPACK_PRL] = 4;
690 return umfpackMatrix_;
699 if (!matrixIsLoaded_)
701 Caller::free_symbolic(&UMF_Symbolic);
702 umfpackMatrix_.free();
704 Caller::free_numeric(&UMF_Numeric);
705 matrixIsLoaded_ =
false;
708 const char* name() {
return "UMFPACK"; }
711 typedef typename Dune::UMFPackMethodChooser<T> Caller;
713 template<
class Mat,
class X,
class TM,
class TD,
class T1>
714 friend class SeqOverlappingSchwarz;
715 friend struct SeqOverlappingSchwarzAssemblerHelper<
UMFPack<
Matrix>,true>;
720 double UMF_Decomposition_Info[UMFPACK_INFO];
721 Caller::symbolic(
static_cast<SuiteSparse_long
>(umfpackMatrix_.N()),
722 static_cast<SuiteSparse_long
>(umfpackMatrix_.N()),
723 umfpackMatrix_.getColStart(),
724 umfpackMatrix_.getRowIndex(),
725 reinterpret_cast<double*
>(umfpackMatrix_.getValues()),
728 UMF_Decomposition_Info);
729 Caller::numeric(umfpackMatrix_.getColStart(),
730 umfpackMatrix_.getRowIndex(),
731 reinterpret_cast<double*
>(umfpackMatrix_.getValues()),
735 UMF_Decomposition_Info);
736 Caller::report_status(UMF_Control,UMF_Decomposition_Info[UMFPACK_STATUS]);
739 std::cout <<
"[UMFPack Decomposition]" << std::endl;
740 std::cout <<
"Wallclock Time taken: " << UMF_Decomposition_Info[UMFPACK_NUMERIC_WALLTIME] <<
" (CPU Time: " << UMF_Decomposition_Info[UMFPACK_NUMERIC_TIME] <<
")" << std::endl;
741 std::cout <<
"Flops taken: " << UMF_Decomposition_Info[UMFPACK_FLOPS] << std::endl;
742 std::cout <<
"Peak Memory Usage: " << UMF_Decomposition_Info[UMFPACK_PEAK_MEMORY]*UMF_Decomposition_Info[UMFPACK_SIZE_OF_UNIT] <<
" bytes" << std::endl;
743 std::cout <<
"Condition number estimate: " << 1./UMF_Decomposition_Info[UMFPACK_RCOND] << std::endl;
744 std::cout <<
"Numbers of non-zeroes in decomposition: L: " << UMF_Decomposition_Info[UMFPACK_LNZ] <<
" U: " << UMF_Decomposition_Info[UMFPACK_UNZ] << std::endl;
748 Caller::report_info(UMF_Control,UMF_Decomposition_Info);
752 void printOnApply(
double* UMF_Info)
754 Caller::report_status(UMF_Control,UMF_Info[UMFPACK_STATUS]);
757 std::cout <<
"[UMFPack Solve]" << std::endl;
758 std::cout <<
"Wallclock Time: " << UMF_Info[UMFPACK_SOLVE_WALLTIME] <<
" (CPU Time: " << UMF_Info[UMFPACK_SOLVE_TIME] <<
")" << std::endl;
759 std::cout <<
"Flops Taken: " << UMF_Info[UMFPACK_SOLVE_FLOPS] << std::endl;
760 std::cout <<
"Iterative Refinement steps taken: " << UMF_Info[UMFPACK_IR_TAKEN] << std::endl;
761 std::cout <<
"Error Estimate: " << UMF_Info[UMFPACK_OMEGA1] <<
" resp. " << UMF_Info[UMFPACK_OMEGA2] << std::endl;
766 bool matrixIsLoaded_;
770 double UMF_Control[UMFPACK_CONTROL];
773 template<
typename T,
typename A,
int n,
int m>
774 struct IsDirectSolver<UMFPack<BCRSMatrix<FieldMatrix<T,n,m>,A> > >
779 template<
typename T,
typename A>
780 struct StoresColumnCompressed<UMFPack<BCRSMatrix<T,A> > >
782 enum { value =
true };
785 struct UMFPackCreator {
786 template<
class F,
class=
void>
struct isValidBlock : std::false_type{};
787 template<
class B>
struct isValidBlock<B, std::enable_if_t<std::is_same<typename FieldTraits<B>::real_type,double>::value>> : std::true_type {};
789 template<
typename TL,
typename M>
790 std::shared_ptr<Dune::InverseOperator<typename Dune::TypeListElement<1, TL>::type,
791 typename Dune::TypeListElement<2, TL>::type>>
794 isValidBlock<
typename Dune::TypeListElement<1, TL>::type::block_type>::value,
int> = 0)
const
796 int verbose = config.
get(
"verbose", 0);
797 return std::make_shared<Dune::UMFPack<M>>(mat,verbose);
801 template<
typename TL,
typename M>
802 std::shared_ptr<Dune::InverseOperator<typename Dune::TypeListElement<1, TL>::type,
803 typename Dune::TypeListElement<2, TL>::type>>
806 !isValidBlock<
typename Dune::TypeListElement<1, TL>::type::block_type>::value,
int> = 0)
const
809 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
812 DUNE_REGISTER_DIRECT_SOLVER(
"umfpack",Dune::UMFPackCreator());
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
Abstract base class for all solvers.
Definition: solver.hh:99
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:185
Default exception class for range errors.
Definition: exceptions.hh:254
The UMFPack direct sparse solver.
Definition: umfpack.hh:245
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_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:697
virtual void apply(domain_type &x, range_type &b, [[maybe_unused]] double reduction, InverseOperatorResult &res)
apply inverse operator, with given convergence criteria.
Definition: umfpack.hh:442
virtual void apply(domain_type &x, range_type &b, InverseOperatorResult &res)
Apply inverse operator,.
Definition: umfpack.hh:391
void * getFactorization()
Return the matrix factorization.
Definition: umfpack.hh:679
virtual SolverCategory::Category category() const
Category of the solver (see SolverCategory::Category)
Definition: umfpack.hh:265
UMFPack(const Matrix &mat_, const ParameterTree &config)
Construct a solver object from a matrix.
Definition: umfpack.hh:315
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:339
typename Impl::UMFPackVectorChooser< M >::range_type range_type
The type of the range of the solver.
Definition: umfpack.hh:262
UMFPack()
default constructor
Definition: umfpack.hh:321
UMFPackMatrix & getInternalMatrix()
Return the column compress matrix from UMFPack.
Definition: umfpack.hh:688
void apply(T *x, T *b)
additional apply method with c-arrays in analogy to superlu
Definition: umfpack.hh:454
void setVerbosity(int v)
sets the verbosity level for the UMFPack solver
Definition: umfpack.hh:663
UMFPack(const char *file, int verbose=0)
try loading a decomposition from file
Definition: umfpack.hh:366
void setMatrix(const Matrix &matrix, const BitVector &bitVector={})
Initialize data from given matrix.
Definition: umfpack.hh:508
void saveDecomposition(const char *file)
saves a decomposition to a file
Definition: umfpack.hh:491
typename Impl::UMFPackVectorChooser< M >::domain_type domain_type
The type of the domain of the solver.
Definition: umfpack.hh:260
UMFPack(const Matrix &matrix, int verbose=0)
Construct a solver object from a matrix.
Definition: umfpack.hh:278
ISTL::Impl::BCCSMatrixInitializer< M, size_type > MatrixInitializer
Type of an associated initializer class.
Definition: umfpack.hh:258
ISTL::Impl::BCCSMatrix< typename Matrix::field_type, size_type > UMFPackMatrix
The corresponding (scalar) UMFPack matrix type.
Definition: umfpack.hh:256
void setOption(unsigned int option, double value)
Set UMFPack-specific options.
Definition: umfpack.hh:480
M Matrix
The matrix type.
Definition: umfpack.hh:253
UMFPack(const Matrix &matrix, int verbose, bool)
Constructor for compatibility with SuperLU standard constructor.
Definition: umfpack.hh:296
Dune namespace.
Definition: alignedallocator.hh:13
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
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
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:48
double elapsed
Elapsed time in seconds.
Definition: solver.hh:82
int iterations
Number of iterations.
Definition: solver.hh:67
bool converged
True if convergence criterion has been met.
Definition: solver.hh:73
Category
Definition: solvercategory.hh:23