3#ifndef DUNE_ISTL_SPQR_HH
4#define DUNE_ISTL_SPQR_HH
6#if HAVE_SUITESPARSE_SPQR || defined DOXYGEN
11#include <SuiteSparseQR.hpp>
15#include <dune/istl/bccsmatrixinitializer.hh>
18#include <dune/istl/solverfactory.hh>
33 template<
class M,
class T,
class TM,
class TD,
class TA>
34 class SeqOverlappingSchwarz;
36 template<
class T,
bool tag>
37 struct SeqOverlappingSchwarzAssemblerHelper;
44 template<
class Matrix>
61 template<
typename T,
typename A,
int n,
int m>
63 :
public InverseOperator<BlockVector<FieldVector<T,m>, typename std::allocator_traits<A>::template rebind_alloc<FieldVector<T,m> > >,
64 BlockVector<FieldVector<T,n>, typename std::allocator_traits<A>::template rebind_alloc<FieldVector<T,n> > > >
73 typedef ISTL::Impl::BCCSMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A>,
int>
MatrixInitializer;
82 return SolverCategory::Category::sequential;
93 SPQR(
const Matrix& matrix,
int verbose=0) : matrixIsLoaded_(false), verbose_(verbose)
96 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
97 "Unsupported Type in SPQR (only double and std::complex<double> supported)");
98 cc_ =
new cholmod_common();
111 SPQR(
const Matrix& matrix,
int verbose,
bool) : matrixIsLoaded_(false), verbose_(verbose)
114 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
115 "Unsupported Type in SPQR (only double and std::complex<double> supported)");
116 cc_ =
new cholmod_common();
117 cholmod_l_start(cc_);
131 :
SPQR(matrix, config.get<int>(
"verbose", 0))
135 SPQR() : matrixIsLoaded_(false), verbose_(0)
138 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
139 "Unsupported Type in SPQR (only double and std::complex<double> supported)");
140 cc_ =
new cholmod_common();
141 cholmod_l_start(cc_);
147 if ((spqrMatrix_.N() + spqrMatrix_.M() > 0) || matrixIsLoaded_)
149 cholmod_l_finish(cc_);
155 const std::size_t dimMat(spqrMatrix_.N());
157 for(std::size_t k = 0; k != dimMat; ++k)
158 (
static_cast<T*
>(B_->x))[k] = b[k];
159 cholmod_dense* BTemp = B_;
160 B_ = SuiteSparseQR_qmult<T>(0, spqrfactorization_, B_, cc_);
161 cholmod_dense* X = SuiteSparseQR_solve<T>(1, spqrfactorization_, B_, cc_);
162 cholmod_l_free_dense(&BTemp, cc_);
164 for(std::size_t k = 0; k != dimMat; ++k)
165 x [k] = (
static_cast<T*
>(X->x))[k];
166 cholmod_l_free_dense(&X, cc_);
172 std::cout<<std::endl<<
"Solving with SuiteSparseQR"<<std::endl;
173 std::cout<<
"Flops Taken: "<<cc_->SPQR_flopcount<<std::endl;
174 std::cout<<
"Analysis Time: "<<cc_->SPQR_analyze_time<<
" s"<<std::endl;
175 std::cout<<
"Factorize Time: "<<cc_->SPQR_factorize_time<<
" s"<<std::endl;
176 std::cout<<
"Backsolve Time: "<<cc_->SPQR_solve_time<<
" s"<<std::endl;
177 std::cout<<
"Peak Memory Usage: "<<cc_->memory_usage<<
" bytes"<<std::endl;
178 std::cout<<
"Rank Estimate: "<<cc_->SPQR_istat[4]<<std::endl<<std::endl;
188 void setOption([[maybe_unused]]
unsigned int option, [[maybe_unused]]
double value)
194 if ((spqrMatrix_.N() + spqrMatrix_.M() > 0) || matrixIsLoaded_)
197 if (spqrMatrix_.N() + spqrMatrix_.M() + spqrMatrix_.nonzeroes() != 0)
199 spqrMatrix_.setSize(MatrixDimension<Matrix>::rowdim(matrix),
200 MatrixDimension<Matrix>::coldim(matrix));
201 ISTL::Impl::BCCSMatrixInitializer<Matrix, int> initializer(spqrMatrix_);
203 copyToBCCSMatrix(initializer, matrix);
209 void setSubMatrix(
const Matrix& matrix,
const S& rowIndexSet)
211 if ((spqrMatrix_.N() + spqrMatrix_.M() > 0) || matrixIsLoaded_)
214 if (spqrMatrix_.N() + spqrMatrix_.M() + spqrMatrix_.nonzeroes() != 0)
217 spqrMatrix_.setSize(rowIndexSet.size()*MatrixDimension<Matrix>::rowdim(matrix) / matrix.
N(),
218 rowIndexSet.size()*MatrixDimension<Matrix>::coldim(matrix) / matrix.
M());
219 ISTL::Impl::BCCSMatrixInitializer<Matrix, int> initializer(spqrMatrix_);
221 copyToBCCSMatrix(initializer, ISTL::Impl::MatrixRowSubset<
Matrix,std::set<std::size_t> >(matrix,rowIndexSet));
241 return spqrfactorization_;
259 cholmod_l_free_sparse(&A_, cc_);
260 cholmod_l_free_dense(&B_, cc_);
261 SuiteSparseQR_free<T>(&spqrfactorization_, cc_);
263 matrixIsLoaded_ =
false;
273 template<
class M,
class X,
class TM,
class TD,
class T1>
276 friend struct SeqOverlappingSchwarzAssemblerHelper<
SPQR<
Matrix>,true>;
281 const std::size_t dimMat(spqrMatrix_.N());
282 const std::size_t nnz(spqrMatrix_.getColStart()[dimMat]);
284 A_ = cholmod_l_allocate_sparse(dimMat, dimMat, nnz, 1, 1, 0, 1, cc_);
286 for(std::size_t k = 0; k != (dimMat+1); ++k)
287 (
static_cast<long int *
>(A_->p))[k] = spqrMatrix_.getColStart()[k];
288 for(std::size_t k = 0; k != nnz; ++k)
290 (
static_cast<long int*
>(A_->i))[k] = spqrMatrix_.getRowIndex()[k];
291 (
static_cast<T*
>(A_->x))[k] = spqrMatrix_.getValues()[k];
294 B_ = cholmod_l_allocate_dense(dimMat, 1, dimMat, A_->xtype, cc_);
296 spqrfactorization_=SuiteSparseQR_factorize<T>(SPQR_ORDERING_DEFAULT,SPQR_DEFAULT_TOL,A_,cc_);
299 SPQRMatrix spqrMatrix_;
300 bool matrixIsLoaded_;
305 SuiteSparseQR_factorization<T>* spqrfactorization_;
308 template<
typename T,
typename A>
309 struct IsDirectSolver<SPQR<BCRSMatrix<T,A> > >
314 template<
typename T,
typename A>
315 struct StoresColumnCompressed<SPQR<BCRSMatrix<T,A> > >
321 template<
class>
struct isValidBlock : std::false_type{};
323 template<
typename TL,
typename M>
324 std::shared_ptr<Dune::InverseOperator<typename Dune::TypeListElement<1, TL>::type,
325 typename Dune::TypeListElement<2, TL>::type>>
328 isValidBlock<
typename Dune::TypeListElement<1, TL>::type::block_type>::value,
int> = 0)
const
330 int verbose = config.
get(
"verbose", 0);
331 return std::make_shared<Dune::SPQR<M>>(mat,verbose);
335 template<
typename TL,
typename M>
336 std::shared_ptr<Dune::InverseOperator<typename Dune::TypeListElement<1, TL>::type,
337 typename Dune::TypeListElement<2, TL>::type>>
339 std::enable_if_t<!isValidBlock<
typename Dune::TypeListElement<1, TL>::type::block_type>::value,
int> = 0)
const
342 "Unsupported Type in SPQR (only double and std::complex<double> supported)");
345 template<>
struct SPQRCreator::isValidBlock<FieldVector<double,1>> : std::true_type{};
348 DUNE_REGISTER_DIRECT_SOLVER(
"spqr", Dune::SPQRCreator());
A sparse block matrix with compressed row storage.
Definition: bcrsmatrix.hh:464
A vector of blocks with memory management.
Definition: bvector.hh:393
A dense n x m matrix.
Definition: fmatrix.hh:69
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
Use the SPQR package to directly solve linear systems – empty default class.
Definition: spqr.hh:46
Sequential overlapping Schwarz preconditioner.
Definition: overlappingschwarz.hh:753
A few common exception classes.
#define DUNE_THROW(E, m)
Definition: exceptions.hh:216
virtual ~SPQR()
Destructor.
Definition: spqr.hh:145
SPQR(const Matrix &matrix, int verbose, bool)
Constructor for compatibility with SuperLU standard constructor.
Definition: spqr.hh:111
virtual SolverCategory::Category category() const
Category of the solver (see SolverCategory::Category)
Definition: spqr.hh:80
SPQRMatrix & getInternalMatrix()
Return the column coppressed matrix.
Definition: spqr.hh:248
void setMatrix(const Matrix &matrix)
Initialize data from given matrix.
Definition: spqr.hh:192
ISTL::Impl::BCCSMatrixInitializer< BCRSMatrix< FieldMatrix< T, n, m >, A >, int > MatrixInitializer
Type of an associated initializer class.
Definition: spqr.hh:73
SPQR()
Default constructor.
Definition: spqr.hh:135
const char * name()
Get method name.
Definition: spqr.hh:267
SuiteSparseQR_factorization< T > * getFactorization()
Return the matrix factorization.
Definition: spqr.hh:239
void setVerbosity(int v)
Sets the verbosity level for the solver.
Definition: spqr.hh:230
ISTL::Impl::BCCSMatrix< T, int > SPQRMatrix
The corresponding SuperLU Matrix type.
Definition: spqr.hh:71
virtual void apply(domain_type &x, range_type &b, double reduction, InverseOperatorResult &res)
apply inverse operator, with given convergence criteria.
Definition: spqr.hh:183
void free()
Free allocated space.
Definition: spqr.hh:257
SPQR(const Matrix &matrix, const ParameterTree &config)
Constructs the SPQR solver.
Definition: spqr.hh:130
Dune::BlockVector< FieldVector< T, n >, typename std::allocator_traits< A >::template rebind_alloc< FieldVector< T, n > > > range_type
The type of the range of the solver.
Definition: spqr.hh:77
Dune::BlockVector< FieldVector< T, m >, typename std::allocator_traits< A >::template rebind_alloc< FieldVector< T, m > > > domain_type
The type of the domain of the solver.
Definition: spqr.hh:75
virtual void apply(domain_type &x, range_type &b, InverseOperatorResult &res)
Apply inverse operator,.
Definition: spqr.hh:153
SPQR(const Matrix &matrix, int verbose=0)
Construct a solver object from a BCRSMatrix.
Definition: spqr.hh:93
Dune namespace.
Definition: alignedallocator.hh:11
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
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