Dune Core Modules (2.6.0)

ldl.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_LDL_HH
4#define DUNE_ISTL_LDL_HH
5
6#if HAVE_SUITESPARSE_LDL || defined DOXYGEN
7
8#include <iostream>
9#include <type_traits>
10
11#ifdef __cplusplus
12extern "C"
13{
14#include "ldl.h"
15#include "amd.h"
16}
17#endif
18
20#include <dune/common/unused.hh>
21
22#include <dune/istl/colcompmatrix.hh>
23#include <dune/istl/solvers.hh>
25
26namespace Dune {
38 // forward declarations
39 template<class M, class T, class TM, class TD, class TA>
40 class SeqOverlappingSchwarz;
41
42 template<class T, bool tag>
43 struct SeqOverlappingSchwarzAssemblerHelper;
44
51 template<class Matrix>
52 class LDL
53 {};
54
68 template<typename T, typename A, int n, int m>
69 class LDL<BCRSMatrix<FieldMatrix<T,n,m>,A > >
70 : public InverseOperator<BlockVector<FieldVector<T,m>, typename A::template rebind<FieldVector<T,m> >::other>,
71 BlockVector<FieldVector<T,n>, typename A::template rebind<FieldVector<T,n> >::other> >
72 {
73 public:
82 typedef Dune::BlockVector<FieldVector<T,m>, typename A::template rebind<FieldVector<T,m> >::other> domain_type;
84 typedef Dune::BlockVector<FieldVector<T,n>, typename A::template rebind<FieldVector<T,n> >::other> range_type;
85
88 {
89 return SolverCategory::Category::sequential;
90 }
91
101 LDL(const Matrix& matrix, int verbose=0) : matrixIsLoaded_(false), verbose_(verbose)
102 {
103 //check whether T is a supported type
104 static_assert(std::is_same<T,double>::value,"Unsupported Type in LDL (only double supported)");
105 setMatrix(matrix);
106 }
107
117 LDL(const Matrix& matrix, int verbose, bool) : matrixIsLoaded_(false), verbose_(verbose)
118 {
119 //check whether T is a supported type
120 static_assert(std::is_same<T,double>::value,"Unsupported Type in LDL (only double supported)");
121 setMatrix(matrix);
122 }
123
125 LDL() : matrixIsLoaded_(false), verbose_(0)
126 {}
127
129 virtual ~LDL()
130 {
131 if ((ldlMatrix_.N() + ldlMatrix_.M() > 0) || matrixIsLoaded_)
132 free();
133 }
134
137 {
138 const int dimMat(ldlMatrix_.N());
139 ldl_perm(dimMat, Y_, reinterpret_cast<double*>(&b[0]), P_);
140 ldl_lsolve(dimMat, Y_, Lp_, Li_, Lx_);
141 ldl_dsolve(dimMat, Y_, D_);
142 ldl_ltsolve(dimMat, Y_, Lp_, Li_, Lx_);
143 ldl_permt(dimMat, reinterpret_cast<double*>(&x[0]), Y_, P_);
144 // this is a direct solver
145 res.iterations = 1;
146 res.converged = true;
147 }
148
150 virtual void apply(domain_type& x, range_type& b, double reduction, InverseOperatorResult& res)
151 {
152 DUNE_UNUSED_PARAMETER(reduction);
153 apply(x,b,res);
154 }
155
161 void apply(T* x, T* b)
162 {
163 const int dimMat(ldlMatrix_.N());
164 ldl_perm(dimMat, Y_, b, P_);
165 ldl_lsolve(dimMat, Y_, Lp_, Li_, Lx_);
166 ldl_dsolve(dimMat, Y_, D_);
167 ldl_ltsolve(dimMat, Y_, Lp_, Li_, Lx_);
168 ldl_permt(dimMat, x, Y_, P_);
169 }
170
171 void setOption(unsigned int option, double value)
172 {
173 DUNE_UNUSED_PARAMETER(option);
175 }
176
178 void setMatrix(const Matrix& matrix)
179 {
180 if ((ldlMatrix_.N() + ldlMatrix_.M() > 0) || matrixIsLoaded_)
181 free();
182 ldlMatrix_ = matrix;
183 decompose();
184 }
185
186 template<class S>
187 void setSubMatrix(const Matrix& matrix, const S& rowIndexSet)
188 {
189 if ((ldlMatrix_.N() + ldlMatrix_.M() > 0) || matrixIsLoaded_)
190 free();
191 ldlMatrix_.setMatrix(matrix,rowIndexSet);
192 decompose();
193 }
194
199 inline void setVerbosity(int v)
200 {
201 verbose_=v;
202 }
203
209 {
210 return ldlMatrix_;
211 }
212
217 void free()
218 {
219 delete [] D_;
220 delete [] Y_;
221 delete [] Lp_;
222 delete [] Lx_;
223 delete [] Li_;
224 delete [] P_;
225 delete [] Pinv_;
226 ldlMatrix_.free();
227 matrixIsLoaded_ = false;
228 }
229
231 inline const char* name()
232 {
233 return "LDL";
234 }
235
240 inline double* getD()
241 {
242 return D_;
243 }
244
249 inline int* getLp()
250 {
251 return Lp_;
252 }
253
258 inline int* getLi()
259 {
260 return Li_;
261 }
262
267 inline double* getLx()
268 {
269 return Lx_;
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<LDL<Matrix>,true>;
277
279 void decompose()
280 {
281 // allocate vectors
282 const int dimMat(ldlMatrix_.N());
283 D_ = new double [dimMat];
284 Y_ = new double [dimMat];
285 Lp_ = new int [dimMat + 1];
286 Parent_ = new int [dimMat];
287 Lnz_ = new int [dimMat];
288 Flag_ = new int [dimMat];
289 Pattern_ = new int [dimMat];
290 P_ = new int [dimMat];
291 Pinv_ = new int [dimMat];
292
293 double Info [AMD_INFO];
294 if(amd_order (dimMat, ldlMatrix_.getColStart(), ldlMatrix_.getRowIndex(), P_, (double *) NULL, Info) < AMD_OK)
295 DUNE_THROW(InvalidStateException,"Error: AMD failed!");
296 if(verbose_ > 0)
297 amd_info (Info);
298 // compute the symbolic factorisation
299 ldl_symbolic(dimMat, ldlMatrix_.getColStart(), ldlMatrix_.getRowIndex(), Lp_, Parent_, Lnz_, Flag_, P_, Pinv_);
300 // initialise those entries of additionalVectors_ whose dimension is known only now
301 Lx_ = new double [Lp_[dimMat]];
302 Li_ = new int [Lp_[dimMat]];
303 // compute the numeric factorisation
304 const int rank(ldl_numeric(dimMat, ldlMatrix_.getColStart(), ldlMatrix_.getRowIndex(), ldlMatrix_.getValues(),
305 Lp_, Parent_, Lnz_, Li_, Lx_, D_, Y_, Pattern_, Flag_, P_, Pinv_));
306 // free temporary vectors
307 delete [] Flag_;
308 delete [] Pattern_;
309 delete [] Parent_;
310 delete [] Lnz_;
311
312 if(rank!=dimMat)
313 DUNE_THROW(InvalidStateException,"Error: LDL factorisation failed!");
314 }
315
316 LDLMatrix ldlMatrix_;
317 bool matrixIsLoaded_;
318 int verbose_;
319 int* Lp_;
320 int* Parent_;
321 int* Lnz_;
322 int* Flag_;
323 int* Pattern_;
324 int* P_;
325 int* Pinv_;
326 double* D_;
327 double* Y_;
328 double* Lx_;
329 int* Li_;
330 };
331
332 template<typename T, typename A, int n, int m>
333 struct IsDirectSolver<LDL<BCRSMatrix<FieldMatrix<T,n,m>,A> > >
334 {
335 enum {value = true};
336 };
337
338 template<typename T, typename A, int n, int m>
339 struct StoresColumnCompressed<LDL<BCRSMatrix<FieldMatrix<T,n,m>,A> > >
340 {
341 enum {value = true};
342 };
343
344}
345
346#endif //HAVE_SUITESPARSE_LDL
347#endif //DUNE_ISTL_LDL_HH
A sparse block matrix with compressed row storage.
Definition: bcrsmatrix.hh:423
A vector of blocks with memory management.
Definition: bvector.hh:317
A dense n x m matrix.
Definition: fmatrix.hh:68
Default exception if a function was called while the object is not in a valid state for that function...
Definition: exceptions.hh:279
Abstract base class for all solvers.
Definition: solver.hh:91
Use the LDL package to directly solve linear systems – empty default class.
Definition: ldl.hh:53
A generic dynamic dense matrix.
Definition: matrix.hh:555
Sequential overlapping Schwarz preconditioner.
Definition: overlappingschwarz.hh:742
A few common exception classes.
decltype(auto) apply(F &&f, ArgTuple &&args)
Apply function with arguments given as tuple.
Definition: apply.hh:58
#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
int * getLp()
Get factorization Lp.
Definition: ldl.hh:249
LDLMatrix & getInternalMatrix()
Return the column compress matrix.
Definition: ldl.hh:208
Dune::BlockVector< FieldVector< T, n >, typename A::template rebind< FieldVector< T, n > >::other > range_type
The type of the range of the solver.
Definition: ldl.hh:84
double * getLx()
Get factorization Lx.
Definition: ldl.hh:267
void setVerbosity(int v)
Sets the verbosity level for the solver.
Definition: ldl.hh:199
void apply(T *x, T *b)
Additional apply method with c-arrays in analogy to superlu.
Definition: ldl.hh:161
ColCompMatrixInitializer< BCRSMatrix< FieldMatrix< T, n, m >, A > > MatrixInitializer
Type of an associated initializer class.
Definition: ldl.hh:80
virtual void apply(domain_type &x, range_type &b, double reduction, InverseOperatorResult &res)
apply inverse operator, with given convergence criteria.
Definition: ldl.hh:150
virtual ~LDL()
Default constructor.
Definition: ldl.hh:129
Dune::BlockVector< FieldVector< T, m >, typename A::template rebind< FieldVector< T, m > >::other > domain_type
The type of the domain of the solver.
Definition: ldl.hh:82
void free()
Free allocated space.
Definition: ldl.hh:217
virtual void apply(domain_type &x, range_type &b, InverseOperatorResult &res)
Apply inverse operator,.
Definition: ldl.hh:136
int * getLi()
Get factorization Li.
Definition: ldl.hh:258
void setMatrix(const Matrix &matrix)
Initialize data from given matrix.
Definition: ldl.hh:178
LDL(const Matrix &matrix, int verbose, bool)
Constructor for compatibility with SuperLU standard constructor.
Definition: ldl.hh:117
double * getD()
Get factorization diagonal matrix D.
Definition: ldl.hh:240
Dune::ColCompMatrix< Matrix > LDLMatrix
The corresponding SuperLU Matrix type.
Definition: ldl.hh:78
virtual SolverCategory::Category category() const
Category of the solver (see SolverCategory::Category)
Definition: ldl.hh:87
const char * name()
Get method name.
Definition: ldl.hh:231
LDL()
Default constructor.
Definition: ldl.hh:125
LDL(const Matrix &matrix, int verbose=0)
Construct a solver object from a BCRSMatrix.
Definition: ldl.hh:101
Dune namespace.
Definition: alignedallocator.hh:10
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:41
int iterations
Number of iterations.
Definition: solver.hh:59
bool converged
True if convergence criterion has been met.
Definition: solver.hh:65
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)