DUNE-FEM (unstable)

umfpack.hh
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright © DUNE Project contributors, see file LICENSE.md in module root
2// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
3// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
4// vi: set et ts=4 sw=2 sts=2:
5#ifndef DUNE_ISTL_UMFPACK_HH
6#define DUNE_ISTL_UMFPACK_HH
7
8#if HAVE_SUITESPARSE_UMFPACK || defined DOXYGEN
9
10#include<complex>
11#include<type_traits>
12
13#include<umfpack.h>
14
18#include<dune/istl/bccsmatrixinitializer.hh>
20#include<dune/istl/matrix.hh>
21#include<dune/istl/foreach.hh>
22#include<dune/istl/multitypeblockmatrix.hh>
23#include<dune/istl/multitypeblockvector.hh>
26#include<dune/istl/solverfactory.hh>
27
28
29
30namespace Dune {
42 // FORWARD DECLARATIONS
43 template<class M, class T, class TM, class TD, class TA>
44 class SeqOverlappingSchwarz;
45
46 template<class T, bool tag>
47 struct SeqOverlappingSchwarzAssemblerHelper;
48
49 // wrapper class for C-Function Calls in the backend. Choose the right function namespace
50 // depending on the template parameter used.
51 template<typename T>
52 struct UMFPackMethodChooser
53 {
54 static constexpr bool valid = false ;
55 };
56
57 template<>
58 struct UMFPackMethodChooser<double>
59 {
60 static constexpr bool valid = true ;
61
62 template<typename... A>
63 static void defaults(A... args)
64 {
65 umfpack_dl_defaults(args...);
66 }
67 template<typename... A>
68 static void free_numeric(A... args)
69 {
70 umfpack_dl_free_numeric(args...);
71 }
72 template<typename... A>
73 static void free_symbolic(A... args)
74 {
75 umfpack_dl_free_symbolic(args...);
76 }
77 template<typename... A>
78 static int load_numeric(A... args)
79 {
80 return umfpack_dl_load_numeric(args...);
81 }
82 template<typename... A>
83 static void numeric(A... args)
84 {
85 umfpack_dl_numeric(args...);
86 }
87 template<typename... A>
88 static void report_info(A... args)
89 {
90 umfpack_dl_report_info(args...);
91 }
92 template<typename... A>
93 static void report_status(A... args)
94 {
95 umfpack_dl_report_status(args...);
96 }
97 template<typename... A>
98 static int save_numeric(A... args)
99 {
100 return umfpack_dl_save_numeric(args...);
101 }
102 template<typename... A>
103 static void solve(A... args)
104 {
105 umfpack_dl_solve(args...);
106 }
107 template<typename... A>
108 static void symbolic(A... args)
109 {
110 umfpack_dl_symbolic(args...);
111 }
112 };
113
114 template<>
115 struct UMFPackMethodChooser<std::complex<double> >
116 {
117 static constexpr bool valid = true ;
118 using size_type = SuiteSparse_long;
119
120 template<typename... A>
121 static void defaults(A... args)
122 {
123 umfpack_zl_defaults(args...);
124 }
125 template<typename... A>
126 static void free_numeric(A... args)
127 {
128 umfpack_zl_free_numeric(args...);
129 }
130 template<typename... A>
131 static void free_symbolic(A... args)
132 {
133 umfpack_zl_free_symbolic(args...);
134 }
135 template<typename... A>
136 static int load_numeric(A... args)
137 {
138 return umfpack_zl_load_numeric(args...);
139 }
140 template<typename... A>
141 static void numeric(const size_type* cs, const size_type* ri, const double* val, A... args)
142 {
143 umfpack_zl_numeric(cs,ri,val,NULL,args...);
144 }
145 template<typename... A>
146 static void report_info(A... args)
147 {
148 umfpack_zl_report_info(args...);
149 }
150 template<typename... A>
151 static void report_status(A... args)
152 {
153 umfpack_zl_report_status(args...);
154 }
155 template<typename... A>
156 static int save_numeric(A... args)
157 {
158 return umfpack_zl_save_numeric(args...);
159 }
160 template<typename... A>
161 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 {
163 const double* cval = reinterpret_cast<const double*>(val);
164 umfpack_zl_solve(m,cs,ri,cval,NULL,x,NULL,b,NULL,args...);
165 }
166 template<typename... A>
167 static void symbolic(size_type m, size_type n, const size_type* cs, const size_type* ri, const double* val, A... args)
168 {
169 umfpack_zl_symbolic(m,n,cs,ri,val,NULL,args...);
170 }
171 };
172
173 namespace Impl
174 {
175 template<class M, class = void>
176 struct UMFPackVectorChooser;
177
179 template<class M> using UMFPackDomainType = typename UMFPackVectorChooser<M>::domain_type;
180
182 template<class M> using UMFPackRangeType = typename UMFPackVectorChooser<M>::range_type;
183
184 template<class M>
185 struct UMFPackVectorChooser<M,
186 std::enable_if_t<(std::is_same<M,double>::value) || (std::is_same<M,std::complex<double> >::value)>>
187 {
188 using domain_type = M;
189 using range_type = M;
190 };
191
192 template<typename T, int n, int m>
193 struct UMFPackVectorChooser<FieldMatrix<T,n,m>,
194 std::enable_if_t<(std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value)>>
195 {
197 using domain_type = FieldVector<T,m>;
199 using range_type = FieldVector<T,n>;
200 };
201
202 template<typename T, typename A>
203 struct UMFPackVectorChooser<BCRSMatrix<T,A>,
204 std::void_t<UMFPackDomainType<T>, UMFPackRangeType<T>>>
205 {
206 // In case of recursive deduction (e.g., BCRSMatrix<FieldMatrix<...>, Allocator<FieldMatrix<...>>>)
207 // the allocator needs to be converted to the sub-block allocator type too (e.g., Allocator<FieldVector<...>>).
208 // Note that matrix allocator is assumed to be the same as the domain/range type of allocators
210 using domain_type = BlockVector<UMFPackDomainType<T>, typename std::allocator_traits<A>::template rebind_alloc<UMFPackDomainType<T>>>;
212 using range_type = BlockVector<UMFPackRangeType<T>, typename std::allocator_traits<A>::template rebind_alloc<UMFPackRangeType<T>>>;
213 };
214
215 template<typename T, typename A>
216 struct UMFPackVectorChooser<Matrix<T,A>,
217 std::void_t<UMFPackDomainType<T>, UMFPackRangeType<T>>>
218 : public UMFPackVectorChooser<BCRSMatrix<T,A>, std::void_t<UMFPackDomainType<T>, UMFPackRangeType<T>>>
219 {};
220
221 // to make the `UMFPackVectorChooser` work with `MultiTypeBlockMatrix`, we need to add an intermediate step for the rows, which are typically `MultiTypeBlockVector`
222 template<typename FirstBlock, typename... Blocks>
223 struct UMFPackVectorChooser<MultiTypeBlockVector<FirstBlock, Blocks...>,
224 std::void_t<UMFPackDomainType<FirstBlock>, UMFPackRangeType<FirstBlock>, UMFPackDomainType<Blocks>...>>
225 {
227 using domain_type = MultiTypeBlockVector<UMFPackDomainType<FirstBlock>, UMFPackDomainType<Blocks>...>;
229 using range_type = UMFPackRangeType<FirstBlock>;
230 };
231
232 // specialization for `MultiTypeBlockMatrix` with `MultiTypeBlockVector` rows
233 template<typename FirstRow, typename... Rows>
234 struct UMFPackVectorChooser<MultiTypeBlockMatrix<FirstRow, Rows...>,
235 std::void_t<UMFPackDomainType<FirstRow>, UMFPackRangeType<FirstRow>, UMFPackRangeType<Rows>...>>
236 {
238 using domain_type = UMFPackDomainType<FirstRow>;
240 using range_type = MultiTypeBlockVector< UMFPackRangeType<FirstRow>, UMFPackRangeType<Rows>... >;
241 };
242
243 // dummy class to represent no BitVector
244 struct NoBitVector
245 {};
246
247
248 }
249
263 template<typename M>
264 class UMFPack : public InverseOperator<Impl::UMFPackDomainType<M>,Impl::UMFPackRangeType<M>>
265 {
266 using T = typename M::field_type;
267
268 public:
269 using size_type = SuiteSparse_long;
270
272 using Matrix = M;
273 using matrix_type = M;
275 using UMFPackMatrix = ISTL::Impl::BCCSMatrix<typename Matrix::field_type, size_type>;
277 using MatrixInitializer = ISTL::Impl::BCCSMatrixInitializer<M, size_type>;
279 using domain_type = Impl::UMFPackDomainType<M>;
281 using range_type = Impl::UMFPackRangeType<M>;
282
285 {
286 return SolverCategory::Category::sequential;
287 }
288
297 UMFPack(const Matrix& matrix, int verbose=0) : matrixIsLoaded_(false)
298 {
299 //check whether T is a supported type
300 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
301 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
302 Caller::defaults(UMF_Control);
303 setVerbosity(verbose);
304 setMatrix(matrix);
305 }
306
315 UMFPack(const Matrix& matrix, int verbose, bool) : matrixIsLoaded_(false)
316 {
317 //check whether T is a supported type
318 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
319 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
320 Caller::defaults(UMF_Control);
321 setVerbosity(verbose);
322 setMatrix(matrix);
323 }
324
343 UMFPack(const Matrix& mat_, const ParameterTree& config) : matrixIsLoaded_(false)
344 {
345 //check whether T is a supported type
346 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
347 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
348 Caller::defaults(UMF_Control);
349 setVerbosity(config.get<int>("verbose", 0));
350 if (config.hasKey("ordering")) {
351 const auto &ordering = config["ordering"];
352 if (ordering == "cholmod") setOption(UMFPACK_ORDERING, UMFPACK_ORDERING_CHOLMOD);
353 else if (ordering == "default") setOption(UMFPACK_ORDERING, UMFPACK_DEFAULT_ORDERING);
354 else if (ordering == "amd") setOption(UMFPACK_ORDERING, UMFPACK_ORDERING_AMD);
355 else if (ordering == "metis") setOption(UMFPACK_ORDERING, UMFPACK_ORDERING_METIS);
356 else if (ordering == "best") setOption(UMFPACK_ORDERING, UMFPACK_ORDERING_BEST);
357 else if (ordering == "none") setOption(UMFPACK_ORDERING, UMFPACK_ORDERING_NONE);
358 #if defined(UMFPACK_VER) && (UMFPACK_VER >= UMFPACK_VER_CODE(6, 0))
359 else if (ordering == "metis_guard") setOption(UMFPACK_ORDERING, UMFPACK_ORDERING_METIS_GUARD);
360 #endif
361 else DUNE_THROW(Dune::NotImplemented, "The UMFPACK ordering '" << ordering << "' is not implemented");
362 }
363 setMatrix(mat_);
364 }
365
368 UMFPack() : matrixIsLoaded_(false), verbosity_(0)
369 {
370 //check whether T is a supported type
371 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
372 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
373 Caller::defaults(UMF_Control);
374 }
375
386 UMFPack(const Matrix& mat_, const char* file, int verbose=0)
387 {
388 //check whether T is a supported type
389 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
390 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
391 Caller::defaults(UMF_Control);
392 setVerbosity(verbose);
393 int errcode = Caller::load_numeric(&UMF_Numeric, const_cast<char*>(file));
394 if ((errcode == UMFPACK_ERROR_out_of_memory) || (errcode == UMFPACK_ERROR_file_IO))
395 {
396 matrixIsLoaded_ = false;
397 setMatrix(mat_);
398 saveDecomposition(file);
399 }
400 else
401 {
402 matrixIsLoaded_ = true;
403 std::cout << "UMFPack decomposition successfully loaded from " << file << std::endl;
404 }
405 }
406
413 UMFPack(const char* file, int verbose=0)
414 {
415 //check whether T is a supported type
416 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
417 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
418 Caller::defaults(UMF_Control);
419 int errcode = Caller::load_numeric(&UMF_Numeric, const_cast<char*>(file));
420 if (errcode == UMFPACK_ERROR_out_of_memory)
421 DUNE_THROW(Dune::Exception, "ran out of memory while loading UMFPack decomposition");
422 if (errcode == UMFPACK_ERROR_file_IO)
423 DUNE_THROW(Dune::Exception, "IO error while loading UMFPack decomposition");
424 matrixIsLoaded_ = true;
425 std::cout << "UMFPack decomposition successfully loaded from " << file << std::endl;
426 setVerbosity(verbose);
427 }
428
429 virtual ~UMFPack()
430 {
431 if ((umfpackMatrix_.N() + umfpackMatrix_.M() > 0) || matrixIsLoaded_)
432 free();
433 }
434
439 {
440 if (umfpackMatrix_.N() != b.dim())
441 DUNE_THROW(Dune::ISTLError, "Size of right-hand-side vector b does not match the number of matrix rows!");
442 if (umfpackMatrix_.M() != x.dim())
443 DUNE_THROW(Dune::ISTLError, "Size of solution vector x does not match the number of matrix columns!");
444 if (b.size() == 0)
445 return;
446
447 // we have to convert x and b into flat structures
448 // however, this is linear in time
449 std::vector<T> xFlat(x.dim()), bFlat(b.dim());
450
451 flatVectorForEach(x, [&](auto&& entry, auto&& offset)
452 {
453 xFlat[ offset ] = entry;
454 });
455
456 flatVectorForEach(b, [&](auto&& entry, auto&& offset)
457 {
458 bFlat[ offset ] = entry;
459 });
460
461 double UMF_Apply_Info[UMFPACK_INFO];
462 Caller::solve(UMFPACK_A,
463 umfpackMatrix_.getColStart(),
464 umfpackMatrix_.getRowIndex(),
465 umfpackMatrix_.getValues(),
466 reinterpret_cast<double*>(&xFlat[0]),
467 reinterpret_cast<double*>(&bFlat[0]),
468 UMF_Numeric,
469 UMF_Control,
470 UMF_Apply_Info);
471
472 // copy back to blocked vector
473 flatVectorForEach(x, [&](auto&& entry, auto offset)
474 {
475 entry = xFlat[offset];
476 });
477
478 //this is a direct solver
479 res.iterations = 1;
480 res.converged = true;
481 res.elapsed = UMF_Apply_Info[UMFPACK_SOLVE_WALLTIME];
482
483 printOnApply(UMF_Apply_Info);
484 }
485
489 void apply (domain_type& x, range_type& b, [[maybe_unused]] double reduction, InverseOperatorResult& res) override
490 {
491 apply(x,b,res);
492 }
493
501 void apply(T* x, T* b)
502 {
503 double UMF_Apply_Info[UMFPACK_INFO];
504 Caller::solve(UMFPACK_A,
505 umfpackMatrix_.getColStart(),
506 umfpackMatrix_.getRowIndex(),
507 umfpackMatrix_.getValues(),
508 x,
509 b,
510 UMF_Numeric,
511 UMF_Control,
512 UMF_Apply_Info);
513 printOnApply(UMF_Apply_Info);
514 }
515
527 void setOption(unsigned int option, double value)
528 {
529 if (option >= UMFPACK_CONTROL)
530 DUNE_THROW(RangeError, "Requested non-existing UMFPack option");
531
532 UMF_Control[option] = value;
533 }
534
538 void saveDecomposition(const char* file)
539 {
540 int errcode = Caller::save_numeric(UMF_Numeric, const_cast<char*>(file));
541 if (errcode != UMFPACK_OK)
542 DUNE_THROW(Dune::Exception,"IO ERROR while trying to save UMFPack decomposition");
543 }
544
554 template<class BitVector = Impl::NoBitVector>
555 void setMatrix(const Matrix& matrix, const BitVector& bitVector = {})
556 {
557 if ((umfpackMatrix_.N() + umfpackMatrix_.M() > 0) || matrixIsLoaded_)
558 free();
559 if (matrix.N() == 0 or matrix.M() == 0)
560 return;
561
562 if (umfpackMatrix_.N() + umfpackMatrix_.M() + umfpackMatrix_.nonzeroes() != 0)
563 umfpackMatrix_.free();
564
565 constexpr bool useBitVector = not std::is_same_v<BitVector,Impl::NoBitVector>;
566
567 // use a dynamic flat vector for the bitset
568 std::vector<bool> flatBitVector;
569 // and a mapping from the compressed indices
570 std::vector<size_type> subIndices;
571
572 [[maybe_unused]] int numberOfIgnoredDofs = 0;
573 int nonZeros = 0;
574
575 if constexpr ( useBitVector )
576 {
577 auto flatSize = flatVectorForEach(bitVector, [](auto&&, auto&&){});
578 flatBitVector.resize(flatSize);
579
580 flatVectorForEach(bitVector, [&](auto&& entry, auto&& offset)
581 {
582 flatBitVector[ offset ] = entry;
583 if ( entry )
584 {
585 numberOfIgnoredDofs++;
586 }
587 });
588 }
589
590 // compute the flat dimension and the number of nonzeros of the matrix
591 auto [flatRows,flatCols] = flatMatrixForEach( matrix, [&](auto&& /*entry*/, auto&& row, auto&& col){
592 // do not count ignored entries
593 if constexpr ( useBitVector )
594 if ( flatBitVector[row] or flatBitVector[col] )
595 return;
596
597 nonZeros++;
598 });
599
600 if constexpr ( useBitVector )
601 {
602 // use the original flatRows!
603 subIndices.resize(flatRows,std::numeric_limits<std::size_t>::max());
604
605 size_type subIndexCounter = 0;
606 for ( size_type i=0; i<size_type(flatRows); i++ )
607 if ( not flatBitVector[ i ] )
608 subIndices[ i ] = subIndexCounter++;
609
610 // update the original matrix size
611 flatRows -= numberOfIgnoredDofs;
612 flatCols -= numberOfIgnoredDofs;
613 }
614
615
616 umfpackMatrix_.setSize(flatRows,flatCols);
617 umfpackMatrix_.Nnz_ = nonZeros;
618
619 // prepare the arrays
620 umfpackMatrix_.colstart = new size_type[flatCols+1];
621 umfpackMatrix_.rowindex = new size_type[nonZeros];
622 umfpackMatrix_.values = new T[nonZeros];
623
624 for ( size_type i=0; i<size_type(flatCols+1); i++ )
625 {
626 umfpackMatrix_.colstart[i] = 0;
627 }
628
629 // at first, we need to compute the column start indices
630 // therefore, we count all entries in each column and in the end we accumulate everything
631 flatMatrixForEach(matrix, [&](auto&& /*entry*/, auto&& flatRowIndex, auto&& flatColIndex)
632 {
633 // do nothing if entry is excluded
634 if constexpr ( useBitVector )
635 if ( flatBitVector[flatRowIndex] or flatBitVector[flatColIndex] )
636 return;
637
638 // pick compressed or uncompressed index
639 // compiler will hopefully do some constexpr optimization here
640 auto colIdx = useBitVector ? subIndices[flatColIndex] : flatColIndex;
641
642 umfpackMatrix_.colstart[colIdx+1]++;
643 });
644
645 // now accumulate
646 for ( size_type i=0; i<(size_type)flatCols; i++ )
647 {
648 umfpackMatrix_.colstart[i+1] += umfpackMatrix_.colstart[i];
649 }
650
651 // we need a compressed position counter in each column
652 std::vector<size_type> colPosition(flatCols,0);
653
654 // now we can set the entries: the procedure below works with both row- or column major index ordering
655 flatMatrixForEach(matrix, [&](auto&& entry, auto&& flatRowIndex, auto&& flatColIndex)
656 {
657 // do nothing if entry is excluded
658 if constexpr ( useBitVector )
659 if ( flatBitVector[flatRowIndex] or flatBitVector[flatColIndex] )
660 return;
661
662 // pick compressed or uncompressed index
663 // compiler will hopefully do some constexpr optimization here
664 auto rowIdx = useBitVector ? subIndices[flatRowIndex] : flatRowIndex;
665 auto colIdx = useBitVector ? subIndices[flatColIndex] : flatColIndex;
666
667 // the start index of each column is already fixed
668 auto colStart = umfpackMatrix_.colstart[colIdx];
669 // get the current number of picked elements in this column
670 auto colPos = colPosition[colIdx];
671 // assign the corresponding row index and the value of this element
672 umfpackMatrix_.rowindex[ colStart + colPos ] = rowIdx;
673 umfpackMatrix_.values[ colStart + colPos ] = entry;
674 // increase the number of picked elements in this column
675 colPosition[colIdx]++;
676 });
677
678 decompose();
679 }
680
681 // Keep legacy version using a set of scalar indices
682 // The new version using a bitVector type for marking the active matrix indices is
683 // directly given in `setMatrix` with an additional BitVector argument.
684 // The new version is more flexible and allows, e.g., marking single components of a matrix block.
685 void setSubMatrix(const Matrix& _mat, const std::set<typename Matrix::size_type>& rowIndexSet)
686 {
687 if ((umfpackMatrix_.N() + umfpackMatrix_.M() > 0) || matrixIsLoaded_)
688 free();
689
690 if (umfpackMatrix_.N() + umfpackMatrix_.M() + umfpackMatrix_.nonzeroes() != 0)
691 umfpackMatrix_.free();
692
693 umfpackMatrix_.setSize(rowIndexSet.size()*MatrixDimension<Matrix>::rowdim(_mat) / _mat.N(),
694 rowIndexSet.size()*MatrixDimension<Matrix>::coldim(_mat) / _mat.M());
695 ISTL::Impl::BCCSMatrixInitializer<Matrix, SuiteSparse_long> initializer(umfpackMatrix_);
696
697 copyToBCCSMatrix(initializer, ISTL::Impl::MatrixRowSubset<Matrix,std::set<typename Matrix::size_type> >(_mat,rowIndexSet));
698
699 decompose();
700 }
701
709 void setVerbosity(int v)
710 {
711 verbosity_ = v;
712 // set the verbosity level in UMFPack
713 if (verbosity_ == 0)
714 UMF_Control[UMFPACK_PRL] = 1;
715 if (verbosity_ == 1)
716 UMF_Control[UMFPACK_PRL] = 2;
717 if (verbosity_ == 2)
718 UMF_Control[UMFPACK_PRL] = 4;
719 }
720
726 {
727 return UMF_Numeric;
728 }
729
735 {
736 return umfpackMatrix_;
737 }
738
743 void free()
744 {
745 if (!matrixIsLoaded_)
746 {
747 Caller::free_symbolic(&UMF_Symbolic);
748 umfpackMatrix_.free();
749 }
750 Caller::free_numeric(&UMF_Numeric);
751 matrixIsLoaded_ = false;
752 }
753
754 const char* name() { return "UMFPACK"; }
755
756 private:
757 typedef typename Dune::UMFPackMethodChooser<T> Caller;
758
759 template<class Mat,class X, class TM, class TD, class T1>
760 friend class SeqOverlappingSchwarz;
761 friend struct SeqOverlappingSchwarzAssemblerHelper<UMFPack<Matrix>,true>;
762
764 void decompose()
765 {
766 double UMF_Decomposition_Info[UMFPACK_INFO];
767 Caller::symbolic(static_cast<SuiteSparse_long>(umfpackMatrix_.N()),
768 static_cast<SuiteSparse_long>(umfpackMatrix_.N()),
769 umfpackMatrix_.getColStart(),
770 umfpackMatrix_.getRowIndex(),
771 reinterpret_cast<double*>(umfpackMatrix_.getValues()),
772 &UMF_Symbolic,
773 UMF_Control,
774 UMF_Decomposition_Info);
775 Caller::numeric(umfpackMatrix_.getColStart(),
776 umfpackMatrix_.getRowIndex(),
777 reinterpret_cast<double*>(umfpackMatrix_.getValues()),
778 UMF_Symbolic,
779 &UMF_Numeric,
780 UMF_Control,
781 UMF_Decomposition_Info);
782 Caller::report_status(UMF_Control,UMF_Decomposition_Info[UMFPACK_STATUS]);
783 if (verbosity_ == 1)
784 {
785 std::cout << "[UMFPack Decomposition]" << std::endl;
786 std::cout << "Wallclock Time taken: " << UMF_Decomposition_Info[UMFPACK_NUMERIC_WALLTIME] << " (CPU Time: " << UMF_Decomposition_Info[UMFPACK_NUMERIC_TIME] << ")" << std::endl;
787 std::cout << "Flops taken: " << UMF_Decomposition_Info[UMFPACK_FLOPS] << std::endl;
788 std::cout << "Peak Memory Usage: " << UMF_Decomposition_Info[UMFPACK_PEAK_MEMORY]*UMF_Decomposition_Info[UMFPACK_SIZE_OF_UNIT] << " bytes" << std::endl;
789 std::cout << "Condition number estimate: " << 1./UMF_Decomposition_Info[UMFPACK_RCOND] << std::endl;
790 std::cout << "Numbers of non-zeroes in decomposition: L: " << UMF_Decomposition_Info[UMFPACK_LNZ] << " U: " << UMF_Decomposition_Info[UMFPACK_UNZ] << std::endl;
791 }
792 if (verbosity_ == 2)
793 {
794 Caller::report_info(UMF_Control,UMF_Decomposition_Info);
795 }
796 }
797
798 void printOnApply(double* UMF_Info)
799 {
800 Caller::report_status(UMF_Control,UMF_Info[UMFPACK_STATUS]);
801 if (verbosity_ > 0)
802 {
803 std::cout << "[UMFPack Solve]" << std::endl;
804 std::cout << "Wallclock Time: " << UMF_Info[UMFPACK_SOLVE_WALLTIME] << " (CPU Time: " << UMF_Info[UMFPACK_SOLVE_TIME] << ")" << std::endl;
805 std::cout << "Flops Taken: " << UMF_Info[UMFPACK_SOLVE_FLOPS] << std::endl;
806 std::cout << "Iterative Refinement steps taken: " << UMF_Info[UMFPACK_IR_TAKEN] << std::endl;
807 std::cout << "Error Estimate: " << UMF_Info[UMFPACK_OMEGA1] << " resp. " << UMF_Info[UMFPACK_OMEGA2] << std::endl;
808 }
809 }
810
811 UMFPackMatrix umfpackMatrix_;
812 bool matrixIsLoaded_;
813 int verbosity_;
814 void *UMF_Symbolic;
815 void *UMF_Numeric;
816 double UMF_Control[UMFPACK_CONTROL];
817 };
818
819 template<typename M>
820 struct IsDirectSolver<UMFPack<M>>
821 {
822 enum { value=true};
823 };
824
825 template<typename T, typename A>
826 struct StoresColumnCompressed<UMFPack<BCRSMatrix<T,A> > >
827 {
828 enum { value = true };
829 };
830
831 namespace UMFPackImpl {
832 template<class OpTraits, class=void> struct isValidBlock : std::false_type{};
833 template<class OpTraits> struct isValidBlock<OpTraits,
834 std::enable_if_t<
835 std::is_same_v<Impl::UMFPackDomainType<typename OpTraits::matrix_type>, typename OpTraits::domain_type>
836 && std::is_same_v<Impl::UMFPackRangeType<typename OpTraits::matrix_type>, typename OpTraits::range_type>
837 && std::is_same_v<typename FieldTraits<typename OpTraits::domain_type::field_type>::real_type, double>
838 && std::is_same_v<typename FieldTraits<typename OpTraits::range_type::field_type>::real_type, double>
839 >> : std::true_type {};
840 }
841 DUNE_REGISTER_SOLVER("umfpack",
842 [](auto opTraits, const auto& op, const Dune::ParameterTree& config)
843 -> std::shared_ptr<typename decltype(opTraits)::solver_type>
844 {
845 using OpTraits = decltype(opTraits);
846 // works only for sequential operators
847 if constexpr (OpTraits::isParallel){
848 if(opTraits.getCommOrThrow(op).communicator().size() > 1)
849 DUNE_THROW(Dune::InvalidStateException, "UMFPack works only for sequential operators.");
850 }
851 if constexpr (OpTraits::isAssembled){
852 using M = typename OpTraits::matrix_type;
853 // check if UMFPack<M>* is convertible to
854 // InverseOperator*. This checks compatibility of the
855 // domain and range types
856 if constexpr (UMFPackImpl::isValidBlock<OpTraits>::value) {
857 const auto& A = opTraits.getAssembledOpOrThrow(op);
858 const M& mat = A->getmat();
859 return std::make_shared<Dune::UMFPack<M>>(mat, config);
860 }
861 }
862 DUNE_THROW(UnsupportedType,
863 "Unsupported Type in UMFPack (only double and std::complex<double> supported)");
864 return nullptr;
865 });
866} // end namespace Dune
867
868#endif // HAVE_SUITESPARSE_UMFPACK
869
870#endif //DUNE_ISTL_UMFPACK_HH
Implementation of the BCRSMatrix class.
Base class for Dune-Exceptions.
Definition: exceptions.hh:98
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:375
Abstract base class for all solvers.
Definition: solver.hh:101
Default exception for dummy implementations.
Definition: exceptions.hh:357
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
bool hasKey(const std::string &key) const
test for key
Definition: parametertree.cc:51
Default exception class for range errors.
Definition: exceptions.hh:348
The UMFPack direct sparse solver.
Definition: umfpack.hh:265
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,...)
Definition: exceptions.hh:314
constexpr auto max
Function object that returns the greater of the given values.
Definition: hybridutilities.hh:489
void free()
free allocated space.
Definition: umfpack.hh:743
UMFPack(const Matrix &mat_, const ParameterTree &config)
Construct a solver object from a matrix.
Definition: umfpack.hh:343
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:386
Impl::UMFPackRangeType< M > range_type
The type of the range of the solver.
Definition: umfpack.hh:281
UMFPack()
default constructor
Definition: umfpack.hh:368
Impl::UMFPackDomainType< M > domain_type
The type of the domain of the solver.
Definition: umfpack.hh:279
void apply(T *x, T *b)
additional apply method with c-arrays in analogy to superlu
Definition: umfpack.hh:501
SolverCategory::Category category() const override
Category of the solver (see SolverCategory::Category)
Definition: umfpack.hh:284
void setVerbosity(int v)
sets the verbosity level for the UMFPack solver
Definition: umfpack.hh:709
UMFPack(const char *file, int verbose=0)
try loading a decomposition from file
Definition: umfpack.hh:413
void apply(domain_type &x, range_type &b, double reduction, InverseOperatorResult &res) override
apply inverse operator, with given convergence criteria.
Definition: umfpack.hh:489
void setMatrix(const Matrix &matrix, const BitVector &bitVector={})
Initialize data from given matrix.
Definition: umfpack.hh:555
void saveDecomposition(const char *file)
saves a decomposition to a file
Definition: umfpack.hh:538
UMFPackMatrix & getInternalMatrix()
Return the column compress matrix from UMFPack.
Definition: umfpack.hh:734
UMFPack(const Matrix &matrix, int verbose=0)
Construct a solver object from a matrix.
Definition: umfpack.hh:297
ISTL::Impl::BCCSMatrixInitializer< M, size_type > MatrixInitializer
Type of an associated initializer class.
Definition: umfpack.hh:277
ISTL::Impl::BCCSMatrix< typename Matrix::field_type, size_type > UMFPackMatrix
The corresponding (scalar) UMFPack matrix type.
Definition: umfpack.hh:275
void setOption(unsigned int option, double value)
Set UMFPack-specific options.
Definition: umfpack.hh:527
M Matrix
The matrix type.
Definition: umfpack.hh:272
void apply(domain_type &x, range_type &b, InverseOperatorResult &res) override
Apply inverse operator,.
Definition: umfpack.hh:438
UMFPack(const Matrix &matrix, int verbose, bool)
Constructor for compatibility with SuperLU standard constructor.
Definition: umfpack.hh:315
void * getFactorization()
Return the matrix factorization.
Definition: umfpack.hh:725
A dynamic dense block matrix class.
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
STL namespace.
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
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden & Uni Heidelberg  |  generated with Hugo v0.111.3 (Mar 8, 23:37, 2026)