Dune Core Modules (2.9.1)

spqr.hh
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (C) 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_SPQR_HH
6#define DUNE_ISTL_SPQR_HH
7
8#if HAVE_SUITESPARSE_SPQR || defined DOXYGEN
9
10#include <complex>
11#include <type_traits>
12
13#include <SuiteSparseQR.hpp>
14
16
17#include <dune/istl/bccsmatrixinitializer.hh>
18#include <dune/istl/solvers.hh>
20#include <dune/istl/solverfactory.hh>
21
22namespace Dune {
34 // forward declarations
35 template<class M, class T, class TM, class TD, class TA>
36 class SeqOverlappingSchwarz;
37
38 template<class T, bool tag>
39 struct SeqOverlappingSchwarzAssemblerHelper;
40
46 template<class Matrix>
47 class SPQR
48 {};
49
63 template<typename T, typename A, int n, int m>
64 class SPQR<BCRSMatrix<FieldMatrix<T,n,m>,A > >
65 : public InverseOperator<BlockVector<FieldVector<T,m>, typename std::allocator_traits<A>::template rebind_alloc<FieldVector<T,m> > >,
66 BlockVector<FieldVector<T,n>, typename std::allocator_traits<A>::template rebind_alloc<FieldVector<T,n> > > >
67 {
68 public:
73 typedef ISTL::Impl::BCCSMatrix<T,int> SPQRMatrix;
75 typedef ISTL::Impl::BCCSMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A>, int> MatrixInitializer;
77 typedef Dune::BlockVector<FieldVector<T,m>, typename std::allocator_traits<A>::template rebind_alloc<FieldVector<T,m> > > domain_type;
79 typedef Dune::BlockVector<FieldVector<T,n>, typename std::allocator_traits<A>::template rebind_alloc<FieldVector<T,n> > > range_type;
80
83 {
84 return SolverCategory::Category::sequential;
85 }
86
95 SPQR(const Matrix& matrix, int verbose=0) : matrixIsLoaded_(false), verbose_(verbose)
96 {
97 //check whether T is a supported type
98 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
99 "Unsupported Type in SPQR (only double and std::complex<double> supported)");
100 cc_ = new cholmod_common();
101 cholmod_l_start(cc_);
102 setMatrix(matrix);
103 }
104
113 SPQR(const Matrix& matrix, int verbose, bool) : matrixIsLoaded_(false), verbose_(verbose)
114 {
115 //check whether T is a supported type
116 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
117 "Unsupported Type in SPQR (only double and std::complex<double> supported)");
118 cc_ = new cholmod_common();
119 cholmod_l_start(cc_);
120 setMatrix(matrix);
121 }
122
132 SPQR(const Matrix& matrix, const ParameterTree& config)
133 : SPQR(matrix, config.get<int>("verbose", 0))
134 {}
135
137 SPQR() : matrixIsLoaded_(false), verbose_(0)
138 {
139 //check whether T is a supported type
140 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
141 "Unsupported Type in SPQR (only double and std::complex<double> supported)");
142 cc_ = new cholmod_common();
143 cholmod_l_start(cc_);
144 }
145
147 virtual ~SPQR()
148 {
149 if ((spqrMatrix_.N() + spqrMatrix_.M() > 0) || matrixIsLoaded_)
150 free();
151 cholmod_l_finish(cc_);
152 }
153
156 {
157 const std::size_t numRows(spqrMatrix_.N());
158 // fill B
159 for(std::size_t k = 0; k != numRows/n; ++k)
160 for (int l = 0; l < n; ++l)
161 (static_cast<T*>(B_->x))[n*k+l] = b[k][l];
162
163 cholmod_dense* BTemp = B_;
164 B_ = SuiteSparseQR_qmult<T>(0, spqrfactorization_, B_, cc_);
165 cholmod_dense* X = SuiteSparseQR_solve<T>(1, spqrfactorization_, B_, cc_);
166 cholmod_l_free_dense(&BTemp, cc_);
167
168 const std::size_t numCols(spqrMatrix_.M());
169 // fill x
170 for(std::size_t k = 0; k != numCols/m; ++k)
171 for (int l = 0; l < m; ++l)
172 x[k][l] = (static_cast<T*>(X->x))[m*k+l];
173
174 cholmod_l_free_dense(&X, cc_);
175 // this is a direct solver
176 res.iterations = 1;
177 res.converged = true;
178 if(verbose_ > 0)
179 {
180 std::cout<<std::endl<<"Solving with SuiteSparseQR"<<std::endl;
181 std::cout<<"Flops Taken: "<<cc_->SPQR_flopcount<<std::endl;
182 std::cout<<"Analysis Time: "<<cc_->SPQR_analyze_time<<" s"<<std::endl;
183 std::cout<<"Factorize Time: "<<cc_->SPQR_factorize_time<<" s"<<std::endl;
184 std::cout<<"Backsolve Time: "<<cc_->SPQR_solve_time<<" s"<<std::endl;
185 std::cout<<"Peak Memory Usage: "<<cc_->memory_usage<<" bytes"<<std::endl;
186 std::cout<<"Rank Estimate: "<<cc_->SPQR_istat[4]<<std::endl<<std::endl;
187 }
188 }
189
191 virtual void apply (domain_type& x, range_type& b, [[maybe_unused]] double reduction, InverseOperatorResult& res)
192 {
193 apply(x, b, res);
194 }
195
196 void setOption([[maybe_unused]] unsigned int option, [[maybe_unused]] double value)
197 {}
198
200 void setMatrix(const Matrix& matrix)
201 {
202 if ((spqrMatrix_.N() + spqrMatrix_.M() > 0) || matrixIsLoaded_)
203 free();
204
205 if (spqrMatrix_.N() + spqrMatrix_.M() + spqrMatrix_.nonzeroes() != 0)
206 spqrMatrix_.free();
207 spqrMatrix_.setSize(MatrixDimension<Matrix>::rowdim(matrix),
208 MatrixDimension<Matrix>::coldim(matrix));
209 ISTL::Impl::BCCSMatrixInitializer<Matrix, int> initializer(spqrMatrix_);
210
211 copyToBCCSMatrix(initializer, matrix);
212
213 decompose();
214 }
215
216 template<class S>
217 void setSubMatrix(const Matrix& matrix, const S& rowIndexSet)
218 {
219 if ((spqrMatrix_.N() + spqrMatrix_.M() > 0) || matrixIsLoaded_)
220 free();
221
222 if (spqrMatrix_.N() + spqrMatrix_.M() + spqrMatrix_.nonzeroes() != 0)
223 spqrMatrix_.free();
224
225 spqrMatrix_.setSize(rowIndexSet.size()*MatrixDimension<Matrix>::rowdim(matrix) / matrix.N(),
226 rowIndexSet.size()*MatrixDimension<Matrix>::coldim(matrix) / matrix.M());
227 ISTL::Impl::BCCSMatrixInitializer<Matrix, int> initializer(spqrMatrix_);
228
229 copyToBCCSMatrix(initializer, ISTL::Impl::MatrixRowSubset<Matrix,std::set<std::size_t> >(matrix,rowIndexSet));
230
231 decompose();
232 }
233
238 inline void setVerbosity(int v)
239 {
240 verbose_=v;
241 }
242
247 inline SuiteSparseQR_factorization<T>* getFactorization()
248 {
249 return spqrfactorization_;
250 }
251
257 {
258 return spqrMatrix_;
259 }
260
265 void free()
266 {
267 cholmod_l_free_sparse(&A_, cc_);
268 cholmod_l_free_dense(&B_, cc_);
269 SuiteSparseQR_free<T>(&spqrfactorization_, cc_);
270 spqrMatrix_.free();
271 matrixIsLoaded_ = false;
272 }
273
275 inline const char* name()
276 {
277 return "SPQR";
278 }
279
280 private:
281 template<class M,class X, class TM, class TD, class T1>
282 friend class SeqOverlappingSchwarz;
283
284 friend struct SeqOverlappingSchwarzAssemblerHelper<SPQR<Matrix>,true>;
285
287 void decompose()
288 {
289 const std::size_t nrows(spqrMatrix_.N());
290 const std::size_t ncols(spqrMatrix_.M());
291 const std::size_t nnz(spqrMatrix_.getColStart()[ncols]);
292
293 // initialise the matrix A (sorted, packed, unsymmetric, real entries)
294 A_ = cholmod_l_allocate_sparse(nrows, ncols, nnz, 1, 1, 0, 1, cc_);
295
296 // copy all the entries of Ap, Ai, Ax
297 for(std::size_t k = 0; k != (ncols+1); ++k)
298 (static_cast<long int *>(A_->p))[k] = spqrMatrix_.getColStart()[k];
299
300 for(std::size_t k = 0; k != nnz; ++k)
301 {
302 (static_cast<long int*>(A_->i))[k] = spqrMatrix_.getRowIndex()[k];
303 (static_cast<T*>(A_->x))[k] = spqrMatrix_.getValues()[k];
304 }
305
306 // initialise the vector B
307 B_ = cholmod_l_allocate_dense(nrows, 1, nrows, A_->xtype, cc_);
308 // compute factorization of A
309 spqrfactorization_=SuiteSparseQR_factorize<T>(SPQR_ORDERING_DEFAULT,SPQR_DEFAULT_TOL,A_,cc_);
310 }
311
312 SPQRMatrix spqrMatrix_;
313 bool matrixIsLoaded_;
314 int verbose_;
315 cholmod_common* cc_;
316 cholmod_sparse* A_;
317 cholmod_dense* B_;
318 SuiteSparseQR_factorization<T>* spqrfactorization_;
319 };
320
321 template<typename T, typename A>
322 struct IsDirectSolver<SPQR<BCRSMatrix<T,A> > >
323 {
324 enum {value = true};
325 };
326
327 template<typename T, typename A>
328 struct StoresColumnCompressed<SPQR<BCRSMatrix<T,A> > >
329 {
330 enum {value = true};
331 };
332
333 struct SPQRCreator {
334 template<class> struct isValidBlock : std::false_type{};
335
336 template<typename TL, typename M>
337 std::shared_ptr<Dune::InverseOperator<typename Dune::TypeListElement<1, TL>::type,
338 typename Dune::TypeListElement<2, TL>::type>>
339 operator() (TL /*tl*/, const M& mat, const Dune::ParameterTree& config,
340 std::enable_if_t<
341 isValidBlock<typename Dune::TypeListElement<1, TL>::type::block_type>::value,int> = 0) const
342 {
343 int verbose = config.get("verbose", 0);
344 return std::make_shared<Dune::SPQR<M>>(mat,verbose);
345 }
346
347 // second version with SFINAE to validate the template parameters of SPQR
348 template<typename TL, typename M>
349 std::shared_ptr<Dune::InverseOperator<typename Dune::TypeListElement<1, TL>::type,
350 typename Dune::TypeListElement<2, TL>::type>>
351 operator() (TL /*tl*/, const M& /*mat*/, const Dune::ParameterTree& /*config*/,
352 std::enable_if_t<!isValidBlock<typename Dune::TypeListElement<1, TL>::type::block_type>::value,int> = 0) const
353 {
354 DUNE_THROW(UnsupportedType,
355 "Unsupported Type in SPQR (only double and std::complex<double> supported)");
356 }
357 };
358 template<> struct SPQRCreator::isValidBlock<FieldVector<double,1>> : std::true_type{};
359 // std::complex is temporary disabled, because it fails if libc++ is used
360 //template<> struct SPQRCreator::isValidMatrixBlock<FieldMatrix<std::complex<double>,1,1>> : std::true_type{};
361 DUNE_REGISTER_DIRECT_SOLVER("spqr", Dune::SPQRCreator());
362
363} // end namespace Dune
364
365
366#endif //HAVE_SUITESPARSE_SPQR
367#endif //DUNE_ISTL_SPQR_HH
A sparse block matrix with compressed row storage.
Definition: bcrsmatrix.hh:466
A vector of blocks with memory management.
Definition: bvector.hh:395
A dense n x m matrix.
Definition: fmatrix.hh:117
Abstract base class for all solvers.
Definition: solver.hh:99
A generic dynamic dense matrix.
Definition: matrix.hh:561
size_type M() const
Return the number of columns.
Definition: matrix.hh:700
size_type N() const
Return the number of rows.
Definition: matrix.hh:695
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
Use the SPQR package to directly solve linear systems – empty default class.
Definition: spqr.hh:48
Sequential overlapping Schwarz preconditioner.
Definition: overlappingschwarz.hh:755
A few common exception classes.
#define DUNE_THROW(E, m)
Definition: exceptions.hh:218
virtual ~SPQR()
Destructor.
Definition: spqr.hh:147
SPQR(const Matrix &matrix, int verbose, bool)
Constructor for compatibility with SuperLU standard constructor.
Definition: spqr.hh:113
virtual SolverCategory::Category category() const
Category of the solver (see SolverCategory::Category)
Definition: spqr.hh:82
SPQRMatrix & getInternalMatrix()
Return the column coppressed matrix.
Definition: spqr.hh:256
void setMatrix(const Matrix &matrix)
Initialize data from given matrix.
Definition: spqr.hh:200
ISTL::Impl::BCCSMatrixInitializer< BCRSMatrix< FieldMatrix< T, n, m >, A >, int > MatrixInitializer
Type of an associated initializer class.
Definition: spqr.hh:75
SPQR()
Default constructor.
Definition: spqr.hh:137
const char * name()
Get method name.
Definition: spqr.hh:275
SuiteSparseQR_factorization< T > * getFactorization()
Return the matrix factorization.
Definition: spqr.hh:247
void setVerbosity(int v)
Sets the verbosity level for the solver.
Definition: spqr.hh:238
ISTL::Impl::BCCSMatrix< T, int > SPQRMatrix
The corresponding SuperLU Matrix type.
Definition: spqr.hh:73
virtual void apply(domain_type &x, range_type &b, double reduction, InverseOperatorResult &res)
apply inverse operator, with given convergence criteria.
Definition: spqr.hh:191
void free()
Free allocated space.
Definition: spqr.hh:265
SPQR(const Matrix &matrix, const ParameterTree &config)
Constructs the SPQR solver.
Definition: spqr.hh:132
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:79
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:77
virtual void apply(domain_type &x, range_type &b, InverseOperatorResult &res)
Apply inverse operator,.
Definition: spqr.hh:155
SPQR(const Matrix &matrix, int verbose=0)
Construct a solver object from a BCRSMatrix.
Definition: spqr.hh:95
Dune namespace.
Definition: alignedallocator.hh:13
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
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
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 21, 23:30, 2024)