Dune Core Modules (2.7.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
14#include <dune/common/unused.hh>
15
16#include <dune/istl/colcompmatrix.hh>
17#include <dune/istl/solvers.hh>
19#include <dune/istl/solverfactory.hh>
20
21namespace Dune {
33 // forward declarations
34 template<class M, class T, class TM, class TD, class TA>
35 class SeqOverlappingSchwarz;
36
37 template<class T, bool tag>
38 struct SeqOverlappingSchwarzAssemblerHelper;
39
45 template<class Matrix>
46 class SPQR
47 {};
48
62 template<typename T, typename A, int n, int m>
63 class SPQR<BCRSMatrix<FieldMatrix<T,n,m>,A > >
64 : public InverseOperator<BlockVector<FieldVector<T,m>, typename A::template rebind<FieldVector<T,m> >::other>,
65 BlockVector<FieldVector<T,n>, typename A::template rebind<FieldVector<T,n> >::other> >
66 {
67 public:
76 typedef Dune::BlockVector<FieldVector<T,m>, typename A::template rebind<FieldVector<T,m> >::other> domain_type;
78 typedef Dune::BlockVector<FieldVector<T,n>, typename A::template rebind<FieldVector<T,n> >::other> range_type;
79
82 {
83 return SolverCategory::Category::sequential;
84 }
85
94 SPQR(const Matrix& matrix, int verbose=0) : matrixIsLoaded_(false), verbose_(verbose)
95 {
96 //check whether T is a supported type
97 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
98 "Unsupported Type in SPQR (only double and std::complex<double> supported)");
99 cc_ = new cholmod_common();
100 cholmod_l_start(cc_);
101 setMatrix(matrix);
102 }
103
112 SPQR(const Matrix& matrix, int verbose, bool) : matrixIsLoaded_(false), verbose_(verbose)
113 {
114 //check whether T is a supported type
115 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
116 "Unsupported Type in SPQR (only double and std::complex<double> supported)");
117 cc_ = new cholmod_common();
118 cholmod_l_start(cc_);
119 setMatrix(matrix);
120 }
121
123 SPQR() : matrixIsLoaded_(false), verbose_(0)
124 {
125 //check whether T is a supported type
126 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
127 "Unsupported Type in SPQR (only double and std::complex<double> supported)");
128 cc_ = new cholmod_common();
129 cholmod_l_start(cc_);
130 }
131
133 virtual ~SPQR()
134 {
135 if ((spqrMatrix_.N() + spqrMatrix_.M() > 0) || matrixIsLoaded_)
136 free();
137 cholmod_l_finish(cc_);
138 }
139
142 {
143 const std::size_t dimMat(spqrMatrix_.N());
144 // fill B
145 for(std::size_t k = 0; k != dimMat; ++k)
146 (static_cast<T*>(B_->x))[k] = b[k];
147 cholmod_dense* BTemp = B_;
148 B_ = SuiteSparseQR_qmult<T>(0, spqrfactorization_, B_, cc_);
149 cholmod_dense* X = SuiteSparseQR_solve<T>(1, spqrfactorization_, B_, cc_);
150 cholmod_l_free_dense(&BTemp, cc_);
151 // fill x
152 for(std::size_t k = 0; k != dimMat; ++k)
153 x [k] = (static_cast<T*>(X->x))[k];
154 cholmod_l_free_dense(&X, cc_);
155 // this is a direct solver
156 res.iterations = 1;
157 res.converged = true;
158 if(verbose_ > 0)
159 {
160 std::cout<<std::endl<<"Solving with SuiteSparseQR"<<std::endl;
161 std::cout<<"Flops Taken: "<<cc_->SPQR_flopcount<<std::endl;
162 std::cout<<"Analysis Time: "<<cc_->SPQR_analyze_time<<" s"<<std::endl;
163 std::cout<<"Factorize Time: "<<cc_->SPQR_factorize_time<<" s"<<std::endl;
164 std::cout<<"Backsolve Time: "<<cc_->SPQR_solve_time<<" s"<<std::endl;
165 std::cout<<"Peak Memory Usage: "<<cc_->memory_usage<<" bytes"<<std::endl;
166 std::cout<<"Rank Estimate: "<<cc_->SPQR_istat[4]<<std::endl<<std::endl;
167 }
168 }
169
171 virtual void apply (domain_type& x, range_type& b, double reduction, InverseOperatorResult& res)
172 {
173 DUNE_UNUSED_PARAMETER(reduction);
174 apply(x, b, res);
175 }
176
177 void setOption(unsigned int option, double value)
178 {
179 DUNE_UNUSED_PARAMETER(option);
181 }
182
184 void setMatrix(const Matrix& matrix)
185 {
186 if ((spqrMatrix_.N() + spqrMatrix_.M() > 0) || matrixIsLoaded_)
187 free();
188 spqrMatrix_ = matrix;
189 decompose();
190 }
191
192 template<class S>
193 void setSubMatrix(const Matrix& matrix, const S& rowIndexSet)
194 {
195 if ((spqrMatrix_.N() + spqrMatrix_.M() > 0) || matrixIsLoaded_)
196 free();
197 spqrMatrix_.setMatrix(matrix,rowIndexSet);
198 decompose();
199 }
200
205 inline void setVerbosity(int v)
206 {
207 verbose_=v;
208 }
209
214 inline SuiteSparseQR_factorization<T>* getFactorization()
215 {
216 return spqrfactorization_;
217 }
218
224 {
225 return spqrMatrix_;
226 }
227
232 void free()
233 {
234 cholmod_l_free_sparse(&A_, cc_);
235 cholmod_l_free_dense(&B_, cc_);
236 SuiteSparseQR_free<T>(&spqrfactorization_, cc_);
237 spqrMatrix_.free();
238 matrixIsLoaded_ = false;
239 }
240
242 inline const char* name()
243 {
244 return "SPQR";
245 }
246
247 private:
248 template<class M,class X, class TM, class TD, class T1>
249 friend class SeqOverlappingSchwarz;
250
251 friend struct SeqOverlappingSchwarzAssemblerHelper<SPQR<Matrix>,true>;
252
254 void decompose()
255 {
256 const std::size_t dimMat(spqrMatrix_.N());
257 const std::size_t nnz(spqrMatrix_.getColStart()[dimMat]);
258 // initialise the matrix A (sorted, packed, unsymmetric, real entries)
259 A_ = cholmod_l_allocate_sparse(dimMat, dimMat, nnz, 1, 1, 0, 1, cc_);
260 // copy all the entries of Ap, Ai, Ax
261 for(std::size_t k = 0; k != (dimMat+1); ++k)
262 (static_cast<long int *>(A_->p))[k] = spqrMatrix_.getColStart()[k];
263 for(std::size_t k = 0; k != nnz; ++k)
264 {
265 (static_cast<long int*>(A_->i))[k] = spqrMatrix_.getRowIndex()[k];
266 (static_cast<T*>(A_->x))[k] = spqrMatrix_.getValues()[k];
267 }
268 // initialise the vector B
269 B_ = cholmod_l_allocate_dense(dimMat, 1, dimMat, A_->xtype, cc_);
270 // compute factorization of A
271 spqrfactorization_=SuiteSparseQR_factorize<T>(SPQR_ORDERING_DEFAULT,SPQR_DEFAULT_TOL,A_,cc_);
272 }
273
274 SPQRMatrix spqrMatrix_;
275 bool matrixIsLoaded_;
276 int verbose_;
277 cholmod_common* cc_;
278 cholmod_sparse* A_;
279 cholmod_dense* B_;
280 SuiteSparseQR_factorization<T>* spqrfactorization_;
281 };
282
283 template<typename T, typename A>
284 struct IsDirectSolver<SPQR<BCRSMatrix<T,A> > >
285 {
286 enum {value = true};
287 };
288
289 template<typename T, typename A>
290 struct StoresColumnCompressed<SPQR<BCRSMatrix<T,A> > >
291 {
292 enum {value = true};
293 };
294
295 struct SPQRCreator {
296 template<class> struct isValidBlock : std::false_type{};
297
298 template<typename TL, typename M>
299 std::shared_ptr<Dune::InverseOperator<typename Dune::TypeListElement<1, TL>::type,
300 typename Dune::TypeListElement<2, TL>::type>>
301 operator() (TL /*tl*/, const M& mat, const Dune::ParameterTree& config,
302 std::enable_if_t<
303 isValidBlock<typename Dune::TypeListElement<1, TL>::type::block_type>::value,int> = 0) const
304 {
305 int verbose = config.get("verbose", 0);
306 return std::make_shared<Dune::SPQR<M>>(mat,verbose);
307 };
308
309 // second version with SFINAE to validate the template parameters of SPQR
310 template<typename TL, typename M>
311 std::shared_ptr<Dune::InverseOperator<typename Dune::TypeListElement<1, TL>::type,
312 typename Dune::TypeListElement<2, TL>::type>>
313 operator() (TL /*tl*/, const M& /*mat*/, const Dune::ParameterTree& /*config*/,
314 std::enable_if_t<!isValidBlock<typename Dune::TypeListElement<1, TL>::type::block_type>::value,int> = 0) const
315 {
316 DUNE_THROW(UnsupportedType,
317 "Unsupported Type in SPQR (only double and std::complex<double> supported)");
318 }
319 };
320 template<> struct SPQRCreator::isValidBlock<FieldVector<double,1>> : std::true_type{};
321 // std::complex is temporary disabled, because it fails if libc++ is used
322 //template<> struct SPQRCreator::isValidMatrixBlock<FieldMatrix<std::complex<double>,1,1>> : std::true_type{};
323 DUNE_REGISTER_DIRECT_SOLVER("spqr", Dune::SPQRCreator());
324
325} // end namespace Dune
326
327
328#endif //HAVE_SUITESPARSE_SPQR
329#endif //DUNE_ISTL_SPQR_HH
A sparse block matrix with compressed row storage.
Definition: bcrsmatrix.hh:425
A vector of blocks with memory management.
Definition: bvector.hh:403
Inititializer for the ColCompMatrix as needed by OverlappingSchwarz.
Definition: colcompmatrix.hh:255
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:558
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:47
Sequential overlapping Schwarz preconditioner.
Definition: overlappingschwarz.hh:751
A few common exception classes.
decltype(auto) apply(F &&f, ArgTuple &&args)
Apply function with arguments given as tuple.
Definition: apply.hh:46
#define DUNE_UNUSED_PARAMETER(parm)
A macro to mark intentionally unused function parameters with.
Definition: unused.hh:25
#define DUNE_THROW(E, m)
Definition: exceptions.hh:216
virtual ~SPQR()
Destructor.
Definition: spqr.hh:133
SPQR(const Matrix &matrix, int verbose, bool)
Constructor for compatibility with SuperLU standard constructor.
Definition: spqr.hh:112
virtual SolverCategory::Category category() const
Category of the solver (see SolverCategory::Category)
Definition: spqr.hh:81
SPQRMatrix & getInternalMatrix()
Return the column coppressed matrix.
Definition: spqr.hh:223
void setMatrix(const Matrix &matrix)
Initialize data from given matrix.
Definition: spqr.hh:184
SPQR()
Default constructor.
Definition: spqr.hh:123
const char * name()
Get method name.
Definition: spqr.hh:242
SuiteSparseQR_factorization< T > * getFactorization()
Return the matrix factorization.
Definition: spqr.hh:214
void setVerbosity(int v)
Sets the verbosity level for the solver.
Definition: spqr.hh:205
Dune::ColCompMatrix< Matrix > SPQRMatrix
The corresponding SuperLU Matrix type.
Definition: spqr.hh:72
Dune::BlockVector< FieldVector< T, m >, typename A::template rebind< FieldVector< T, m > >::other > domain_type
The type of the domain of the solver.
Definition: spqr.hh:76
virtual void apply(domain_type &x, range_type &b, double reduction, InverseOperatorResult &res)
apply inverse operator, with given convergence criteria.
Definition: spqr.hh:171
void free()
Free allocated space.
Definition: spqr.hh:232
Dune::BlockVector< FieldVector< T, n >, typename A::template rebind< FieldVector< T, n > >::other > range_type
The type of the range of the solver.
Definition: spqr.hh:78
virtual void apply(domain_type &x, range_type &b, InverseOperatorResult &res)
Apply inverse operator,.
Definition: spqr.hh:141
SPQR(const Matrix &matrix, int verbose=0)
Construct a solver object from a BCRSMatrix.
Definition: spqr.hh:94
Dune namespace.
Definition: alignedallocator.hh:14
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
Definition of the DUNE_UNUSED macro for the case that config.h is not available.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Jul 15, 22:36, 2024)