Dune Core Modules (2.5.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
20namespace 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 A::template rebind<FieldVector<T,m> >::other>,
64 BlockVector<FieldVector<T,n>, typename A::template rebind<FieldVector<T,n> >::other> >
65 {
66 public:
75 typedef Dune::BlockVector<FieldVector<T,m>, typename A::template rebind<FieldVector<T,m> >::other> domain_type;
77 typedef Dune::BlockVector<FieldVector<T,n>, typename A::template rebind<FieldVector<T,n> >::other> range_type;
78
87 SPQR(const Matrix& matrix, int verbose=0) : matrixIsLoaded_(false), verbose_(verbose)
88 {
89 //check whether T is a supported type
90 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
91 "Unsupported Type in SPQR (only double and std::complex<double> supported)");
92 cc_ = new cholmod_common();
93 cholmod_l_start(cc_);
94 setMatrix(matrix);
95 }
96
105 SPQR(const Matrix& matrix, int verbose, bool) : matrixIsLoaded_(false), verbose_(verbose)
106 {
107 //check whether T is a supported type
108 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
109 "Unsupported Type in SPQR (only double and std::complex<double> supported)");
110 cc_ = new cholmod_common();
111 cholmod_l_start(cc_);
112 setMatrix(matrix);
113 }
114
116 SPQR() : matrixIsLoaded_(false), verbose_(0)
117 {
118 //check whether T is a supported type
119 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
120 "Unsupported Type in SPQR (only double and std::complex<double> supported)");
121 cc_ = new cholmod_common();
122 cholmod_l_start(cc_);
123 }
124
126 virtual ~SPQR()
127 {
128 if ((spqrMatrix_.N() + spqrMatrix_.M() > 0) || matrixIsLoaded_)
129 free();
130 cholmod_l_finish(cc_);
131 }
132
135 {
136 const std::size_t dimMat(spqrMatrix_.N());
137 // fill B
138 for(std::size_t k = 0; k != dimMat; ++k)
139 (static_cast<T*>(B_->x))[k] = b[k];
140 cholmod_dense* BTemp = B_;
141 B_ = SuiteSparseQR_qmult<T>(0, spqrfactorization_, B_, cc_);
142 cholmod_dense* X = SuiteSparseQR_solve<T>(1, spqrfactorization_, B_, cc_);
143 cholmod_l_free_dense(&BTemp, cc_);
144 // fill x
145 for(std::size_t k = 0; k != dimMat; ++k)
146 x [k] = (static_cast<T*>(X->x))[k];
147 cholmod_l_free_dense(&X, cc_);
148 // this is a direct solver
149 res.iterations = 1;
150 res.converged = true;
151 if(verbose_ > 0)
152 {
153 std::cout<<std::endl<<"Solving with SuiteSparseQR"<<std::endl;
154 std::cout<<"Flops Taken: "<<cc_->SPQR_flopcount<<std::endl;
155 std::cout<<"Analysis Time: "<<cc_->SPQR_analyze_time<<" s"<<std::endl;
156 std::cout<<"Factorize Time: "<<cc_->SPQR_factorize_time<<" s"<<std::endl;
157 std::cout<<"Backsolve Time: "<<cc_->SPQR_solve_time<<" s"<<std::endl;
158 std::cout<<"Peak Memory Usage: "<<cc_->memory_usage<<" bytes"<<std::endl;
159 std::cout<<"Rank Estimate: "<<cc_->SPQR_istat[4]<<std::endl<<std::endl;
160 }
161 }
162
164 virtual void apply (domain_type& x, range_type& b, double reduction, InverseOperatorResult& res)
165 {
166 DUNE_UNUSED_PARAMETER(reduction);
167 apply(x, b, res);
168 }
169
170 void setOption(unsigned int option, double value)
171 {
172 DUNE_UNUSED_PARAMETER(option);
174 }
175
177 void setMatrix(const Matrix& matrix)
178 {
179 if ((spqrMatrix_.N() + spqrMatrix_.M() > 0) || matrixIsLoaded_)
180 free();
181 spqrMatrix_ = matrix;
182 decompose();
183 }
184
185 template<class S>
186 void setSubMatrix(const Matrix& matrix, const S& rowIndexSet)
187 {
188 if ((spqrMatrix_.N() + spqrMatrix_.M() > 0) || matrixIsLoaded_)
189 free();
190 spqrMatrix_.setMatrix(matrix,rowIndexSet);
191 decompose();
192 }
193
198 inline void setVerbosity(int v)
199 {
200 verbose_=v;
201 }
202
207 inline SuiteSparseQR_factorization<T>* getFactorization()
208 {
209 return spqrfactorization_;
210 }
211
217 {
218 return spqrMatrix_;
219 }
220
225 void free()
226 {
227 cholmod_l_free_sparse(&A_, cc_);
228 cholmod_l_free_dense(&B_, cc_);
229 SuiteSparseQR_free<T>(&spqrfactorization_, cc_);
230 spqrMatrix_.free();
231 matrixIsLoaded_ = false;
232 }
233
235 inline const char* name()
236 {
237 return "SPQR";
238 }
239
240 private:
241 template<class M,class X, class TM, class TD, class T1>
242 friend class SeqOverlappingSchwarz;
243
244 friend struct SeqOverlappingSchwarzAssemblerHelper<SPQR<Matrix>,true>;
245
247 void decompose()
248 {
249 const std::size_t dimMat(spqrMatrix_.N());
250 const std::size_t nnz(spqrMatrix_.getColStart()[dimMat]);
251 // initialise the matrix A (sorted, packed, unsymmetric, real entries)
252 A_ = cholmod_l_allocate_sparse(dimMat, dimMat, nnz, 1, 1, 0, 1, cc_);
253 // copy all the entries of Ap, Ai, Ax
254 for(std::size_t k = 0; k != (dimMat+1); ++k)
255 (static_cast<long int *>(A_->p))[k] = spqrMatrix_.getColStart()[k];
256 for(std::size_t k = 0; k != nnz; ++k)
257 {
258 (static_cast<long int*>(A_->i))[k] = spqrMatrix_.getRowIndex()[k];
259 (static_cast<T*>(A_->x))[k] = spqrMatrix_.getValues()[k];
260 }
261 // initialise the vector B
262 B_ = cholmod_l_allocate_dense(dimMat, 1, dimMat, A_->xtype, cc_);
263 // compute factorization of A
264 spqrfactorization_=SuiteSparseQR_factorize<T>(SPQR_ORDERING_DEFAULT,SPQR_DEFAULT_TOL,A_,cc_);
265 }
266
267 SPQRMatrix spqrMatrix_;
268 bool matrixIsLoaded_;
269 int verbose_;
270 cholmod_common* cc_;
271 cholmod_sparse* A_;
272 cholmod_dense* B_;
273 SuiteSparseQR_factorization<T>* spqrfactorization_;
274 };
275
276 template<typename T, typename A, int n, int m>
277 struct IsDirectSolver<SPQR<BCRSMatrix<FieldMatrix<T,n,m>,A> > >
278 {
279 enum {value = true};
280 };
281
282 template<typename T, typename A, int n, int m>
283 struct StoresColumnCompressed<SPQR<BCRSMatrix<FieldMatrix<T,n,m>,A> > >
284 {
285 enum {value = true};
286 };
287
288}
289
290#endif //HAVE_SUITESPARSE_SPQR
291#endif //DUNE_ISTL_SPQR_HH
A sparse block matrix with compressed row storage.
Definition: bcrsmatrix.hh:423
A vector of blocks with memory management.
Definition: bvector.hh:313
A dense n x m matrix.
Definition: fmatrix.hh:68
Abstract base class for all solvers.
Definition: solver.hh:79
A generic dynamic dense matrix.
Definition: matrix.hh:555
Use the SPQR package to directly solve linear systems – empty default class.
Definition: spqr.hh:46
Sequential overlapping Schwarz preconditioner.
Definition: overlappingschwarz.hh:742
A few common exception classes.
virtual ~SPQR()
Destructor.
Definition: spqr.hh:126
SPQR(const Matrix &matrix, int verbose, bool)
Constructor for compatibility with SuperLU standard constructor.
Definition: spqr.hh:105
SPQRMatrix & getInternalMatrix()
Return the column coppressed matrix.
Definition: spqr.hh:216
void setMatrix(const Matrix &matrix)
Initialize data from given matrix.
Definition: spqr.hh:177
SPQR()
Default constructor.
Definition: spqr.hh:116
const char * name()
Get method name.
Definition: spqr.hh:235
SuiteSparseQR_factorization< T > * getFactorization()
Return the matrix factorization.
Definition: spqr.hh:207
void setVerbosity(int v)
Sets the verbosity level for the solver.
Definition: spqr.hh:198
Dune::ColCompMatrix< Matrix > SPQRMatrix
The corresponding SuperLU Matrix type.
Definition: spqr.hh:71
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:75
virtual void apply(domain_type &x, range_type &b, double reduction, InverseOperatorResult &res)
apply inverse operator, with given convergence criteria.
Definition: spqr.hh:164
void free()
Free allocated space.
Definition: spqr.hh:225
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:77
ColCompMatrixInitializer< BCRSMatrix< FieldMatrix< T, n, m >, A > > MatrixInitializer
Type of an associated initializer class.
Definition: spqr.hh:73
virtual void apply(domain_type &x, range_type &b, InverseOperatorResult &res)
Apply inverse operator,.
Definition: spqr.hh:134
SPQR(const Matrix &matrix, int verbose=0)
Construct a solver object from a BCRSMatrix.
Definition: spqr.hh:87
Dune namespace.
Definition: alignment.hh:11
Implementations of the inverse operator interface.
Templates characterizing the type of a solver.
Inititializer for the ColCompMatrix as needed by OverlappingSchwarz.
Definition: colcompmatrix.hh:154
Statistics about the application of an inverse operator.
Definition: solver.hh:32
int iterations
Number of iterations.
Definition: solver.hh:50
bool converged
True if convergence criterion has been met.
Definition: solver.hh:56
Definition of the DUNE_UNUSED macro for the case that config.h is not available.
#define DUNE_UNUSED_PARAMETER(parm)
A macro to mark intentionally unused function parameters with.
Definition: unused.hh:18
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden  |  generated with Hugo v0.111.3 (Nov 12, 23:30, 2024)