dune-istl  2.8.0
spqr.hh
Go to the documentation of this file.
1 // -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2 // vi: set et ts=4 sw=2 sts=2:
3 #ifndef DUNE_ISTL_SPQR_HH
4 #define DUNE_ISTL_SPQR_HH
5 
6 #if HAVE_SUITESPARSE_SPQR || defined DOXYGEN
7 
8 #include <complex>
9 #include <type_traits>
10 
11 #include <SuiteSparseQR.hpp>
12 
13 #include <dune/common/exceptions.hh>
14 
16 #include <dune/istl/solvers.hh>
17 #include <dune/istl/solvertype.hh>
19 
20 namespace Dune {
32  // forward declarations
33  template<class M, class T, class TM, class TD, class TA>
34  class SeqOverlappingSchwarz;
35 
36  template<class T, bool tag>
37  struct SeqOverlappingSchwarzAssemblerHelper;
38 
44  template<class Matrix>
45  class SPQR
46  {};
47 
61  template<typename T, typename A, int n, int m>
62  class SPQR<BCRSMatrix<FieldMatrix<T,n,m>,A > >
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> > > >
65  {
66  public:
71  typedef ISTL::Impl::BCCSMatrix<T,int> SPQRMatrix;
73  typedef ISTL::Impl::BCCSMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A>, int> MatrixInitializer;
75  typedef Dune::BlockVector<FieldVector<T,m>, typename std::allocator_traits<A>::template rebind_alloc<FieldVector<T,m> > > domain_type;
77  typedef Dune::BlockVector<FieldVector<T,n>, typename std::allocator_traits<A>::template rebind_alloc<FieldVector<T,n> > > range_type;
78 
81  {
82  return SolverCategory::Category::sequential;
83  }
84 
93  SPQR(const Matrix& matrix, int verbose=0) : matrixIsLoaded_(false), verbose_(verbose)
94  {
95  //check whether T is a supported type
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();
99  cholmod_l_start(cc_);
100  setMatrix(matrix);
101  }
102 
111  SPQR(const Matrix& matrix, int verbose, bool) : matrixIsLoaded_(false), verbose_(verbose)
112  {
113  //check whether T is a supported type
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_);
118  setMatrix(matrix);
119  }
120 
130  SPQR(const Matrix& matrix, const ParameterTree& config)
131  : SPQR(matrix, config.get<int>("verbose", 0))
132  {}
133 
135  SPQR() : matrixIsLoaded_(false), verbose_(0)
136  {
137  //check whether T is a supported type
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_);
142  }
143 
145  virtual ~SPQR()
146  {
147  if ((spqrMatrix_.N() + spqrMatrix_.M() > 0) || matrixIsLoaded_)
148  free();
149  cholmod_l_finish(cc_);
150  }
151 
154  {
155  const std::size_t dimMat(spqrMatrix_.N());
156  // fill B
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_);
163  // fill x
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_);
167  // this is a direct solver
168  res.iterations = 1;
169  res.converged = true;
170  if(verbose_ > 0)
171  {
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;
179  }
180  }
181 
183  virtual void apply (domain_type& x, range_type& b, [[maybe_unused]] double reduction, InverseOperatorResult& res)
184  {
185  apply(x, b, res);
186  }
187 
188  void setOption([[maybe_unused]] unsigned int option, [[maybe_unused]] double value)
189  {}
190 
192  void setMatrix(const Matrix& matrix)
193  {
194  if ((spqrMatrix_.N() + spqrMatrix_.M() > 0) || matrixIsLoaded_)
195  free();
196 
197  if (spqrMatrix_.N() + spqrMatrix_.M() + spqrMatrix_.nonzeroes() != 0)
198  spqrMatrix_.free();
199  spqrMatrix_.setSize(MatrixDimension<Matrix>::rowdim(matrix),
201  ISTL::Impl::BCCSMatrixInitializer<Matrix, int> initializer(spqrMatrix_);
202 
203  copyToBCCSMatrix(initializer, matrix);
204 
205  decompose();
206  }
207 
208  template<class S>
209  void setSubMatrix(const Matrix& matrix, const S& rowIndexSet)
210  {
211  if ((spqrMatrix_.N() + spqrMatrix_.M() > 0) || matrixIsLoaded_)
212  free();
213 
214  if (spqrMatrix_.N() + spqrMatrix_.M() + spqrMatrix_.nonzeroes() != 0)
215  spqrMatrix_.free();
216 
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_);
220 
221  copyToBCCSMatrix(initializer, ISTL::Impl::MatrixRowSubset<Matrix,std::set<std::size_t> >(matrix,rowIndexSet));
222 
223  decompose();
224  }
225 
230  inline void setVerbosity(int v)
231  {
232  verbose_=v;
233  }
234 
239  inline SuiteSparseQR_factorization<T>* getFactorization()
240  {
241  return spqrfactorization_;
242  }
243 
249  {
250  return spqrMatrix_;
251  }
252 
257  void free()
258  {
259  cholmod_l_free_sparse(&A_, cc_);
260  cholmod_l_free_dense(&B_, cc_);
261  SuiteSparseQR_free<T>(&spqrfactorization_, cc_);
262  spqrMatrix_.free();
263  matrixIsLoaded_ = false;
264  }
265 
267  inline const char* name()
268  {
269  return "SPQR";
270  }
271 
272  private:
273  template<class M,class X, class TM, class TD, class T1>
274  friend class SeqOverlappingSchwarz;
275 
276  friend struct SeqOverlappingSchwarzAssemblerHelper<SPQR<Matrix>,true>;
277 
279  void decompose()
280  {
281  const std::size_t dimMat(spqrMatrix_.N());
282  const std::size_t nnz(spqrMatrix_.getColStart()[dimMat]);
283  // initialise the matrix A (sorted, packed, unsymmetric, real entries)
284  A_ = cholmod_l_allocate_sparse(dimMat, dimMat, nnz, 1, 1, 0, 1, cc_);
285  // copy all the entries of Ap, Ai, Ax
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)
289  {
290  (static_cast<long int*>(A_->i))[k] = spqrMatrix_.getRowIndex()[k];
291  (static_cast<T*>(A_->x))[k] = spqrMatrix_.getValues()[k];
292  }
293  // initialise the vector B
294  B_ = cholmod_l_allocate_dense(dimMat, 1, dimMat, A_->xtype, cc_);
295  // compute factorization of A
296  spqrfactorization_=SuiteSparseQR_factorize<T>(SPQR_ORDERING_DEFAULT,SPQR_DEFAULT_TOL,A_,cc_);
297  }
298 
299  SPQRMatrix spqrMatrix_;
300  bool matrixIsLoaded_;
301  int verbose_;
302  cholmod_common* cc_;
303  cholmod_sparse* A_;
304  cholmod_dense* B_;
305  SuiteSparseQR_factorization<T>* spqrfactorization_;
306  };
307 
308  template<typename T, typename A>
310  {
311  enum {value = true};
312  };
313 
314  template<typename T, typename A>
316  {
317  enum {value = true};
318  };
319 
320  struct SPQRCreator {
321  template<class> struct isValidBlock : std::false_type{};
322 
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>>
326  operator() (TL /*tl*/, const M& mat, const Dune::ParameterTree& config,
327  std::enable_if_t<
328  isValidBlock<typename Dune::TypeListElement<1, TL>::type::block_type>::value,int> = 0) const
329  {
330  int verbose = config.get("verbose", 0);
331  return std::make_shared<Dune::SPQR<M>>(mat,verbose);
332  }
333 
334  // second version with SFINAE to validate the template parameters of SPQR
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>>
338  operator() (TL /*tl*/, const M& /*mat*/, const Dune::ParameterTree& /*config*/,
339  std::enable_if_t<!isValidBlock<typename Dune::TypeListElement<1, TL>::type::block_type>::value,int> = 0) const
340  {
341  DUNE_THROW(UnsupportedType,
342  "Unsupported Type in SPQR (only double and std::complex<double> supported)");
343  }
344  };
345  template<> struct SPQRCreator::isValidBlock<FieldVector<double,1>> : std::true_type{};
346  // std::complex is temporary disabled, because it fails if libc++ is used
347  //template<> struct SPQRCreator::isValidMatrixBlock<FieldMatrix<std::complex<double>,1,1>> : std::true_type{};
349 
350 } // end namespace Dune
351 
352 
353 #endif //HAVE_SUITESPARSE_SPQR
354 #endif //DUNE_ISTL_SPQR_HH
Implementations of the inverse operator interface.
Templates characterizing the type of a solver.
virtual ~SPQR()
Destructor.
Definition: spqr.hh:145
SPQRMatrix & getInternalMatrix()
Return the column coppressed matrix.
Definition: spqr.hh:248
const char * name()
Get method name.
Definition: spqr.hh:267
void setOption([[maybe_unused]] unsigned int option, [[maybe_unused]] double value)
Definition: spqr.hh:188
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
void setMatrix(const Matrix &matrix)
Initialize data from given matrix.
Definition: spqr.hh:192
SuiteSparseQR_factorization< T > * getFactorization()
Return the matrix factorization.
Definition: spqr.hh:239
DUNE_REGISTER_DIRECT_SOLVER("ldl", Dune::LDLCreator())
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
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
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
virtual void apply(domain_type &x, range_type &b, [[maybe_unused]] double reduction, InverseOperatorResult &res)
apply inverse operator, with given convergence criteria.
Definition: spqr.hh:183
void setSubMatrix(const Matrix &matrix, const S &rowIndexSet)
Definition: spqr.hh:209
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
std::shared_ptr< Dune::InverseOperator< typename Dune::TypeListElement< 1, TL >::type, typename Dune::TypeListElement< 2, TL >::type > > operator()(TL, const M &mat, const Dune::ParameterTree &config, std::enable_if_t< isValidBlock< typename Dune::TypeListElement< 1, TL >::type::block_type >::value, int >=0) const
Definition: spqr.hh:326
Matrix & mat
Definition: matrixmatrix.hh:345
Definition: allocator.hh:9
PropertyMapTypeSelector< Amg::VertexVisitedTag, Amg::PropertiesGraph< G, Amg::VertexProperties, EP, VM, EM > >::Type get([[maybe_unused]] const Amg::VertexVisitedTag &tag, Amg::PropertiesGraph< G, Amg::VertexProperties, EP, VM, EM > &graph)
Definition: dependency.hh:291
Definition: matrixutils.hh:209
A sparse block matrix with compressed row storage.
Definition: bcrsmatrix.hh:464
size_type M() const
number of columns (counted in blocks)
Definition: bcrsmatrix.hh:1976
size_type N() const
number of rows (counted in blocks)
Definition: bcrsmatrix.hh:1970
A vector of blocks with memory management.
Definition: bvector.hh:393
Sequential overlapping Schwarz preconditioner.
Definition: overlappingschwarz.hh:753
Definition: overlappingschwarz.hh:692
Definition: matrixutils.hh:25
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
Abstract base class for all solvers.
Definition: solver.hh:97
Category
Definition: solvercategory.hh:21
Definition: solverregistry.hh:75
Definition: solvertype.hh:14
@ value
Whether this is a direct solver.
Definition: solvertype.hh:22
Definition: solvertype.hh:28
@ value
whether the solver internally uses column compressed storage
Definition: solvertype.hh:34
Use the SPQR package to directly solve linear systems – empty default class.
Definition: spqr.hh:46
Definition: spqr.hh:320
Definition: spqr.hh:321