DUNE PDELab (2.7)

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 <memory>
10#include <type_traits>
11
12#ifdef __cplusplus
13extern "C"
14{
15#include "ldl.h"
16#include "amd.h"
17}
18#endif
19
21#include <dune/common/unused.hh>
22
23#include <dune/istl/colcompmatrix.hh>
24#include <dune/istl/solvers.hh>
26#include <dune/istl/solverfactory.hh>
27
28namespace Dune {
40 // forward declarations
41 template<class M, class T, class TM, class TD, class TA>
42 class SeqOverlappingSchwarz;
43
44 template<class T, bool tag>
45 struct SeqOverlappingSchwarzAssemblerHelper;
46
53 template<class Matrix>
54 class LDL
55 {};
56
70 template<typename T, typename A, int n, int m>
71 class LDL<BCRSMatrix<FieldMatrix<T,n,m>,A > >
72 : public InverseOperator<BlockVector<FieldVector<T,m>, typename std::allocator_traits<A>::template rebind_alloc<FieldVector<T,m> > >,
73 BlockVector<FieldVector<T,n>, typename std::allocator_traits<A>::template rebind_alloc<FieldVector<T,n> > > >
74 {
75 public:
84 typedef Dune::BlockVector<FieldVector<T,m>, typename std::allocator_traits<A>::template rebind_alloc<FieldVector<T,m> > > domain_type;
86 typedef Dune::BlockVector<FieldVector<T,n>, typename std::allocator_traits<A>::template rebind_alloc<FieldVector<T,n> > > range_type;
87
90 {
91 return SolverCategory::Category::sequential;
92 }
93
103 LDL(const Matrix& matrix, int verbose=0) : matrixIsLoaded_(false), verbose_(verbose)
104 {
105 //check whether T is a supported type
106 static_assert(std::is_same<T,double>::value,"Unsupported Type in LDL (only double supported)");
107 setMatrix(matrix);
108 }
109
119 LDL(const Matrix& matrix, int verbose, bool) : matrixIsLoaded_(false), verbose_(verbose)
120 {
121 //check whether T is a supported type
122 static_assert(std::is_same<T,double>::value,"Unsupported Type in LDL (only double supported)");
123 setMatrix(matrix);
124 }
125
127 LDL() : matrixIsLoaded_(false), verbose_(0)
128 {}
129
131 virtual ~LDL()
132 {
133 if ((ldlMatrix_.N() + ldlMatrix_.M() > 0) || matrixIsLoaded_)
134 free();
135 }
136
139 {
140 const int dimMat(ldlMatrix_.N());
141 ldl_perm(dimMat, Y_, reinterpret_cast<double*>(&b[0]), P_);
142 ldl_lsolve(dimMat, Y_, Lp_, Li_, Lx_);
143 ldl_dsolve(dimMat, Y_, D_);
144 ldl_ltsolve(dimMat, Y_, Lp_, Li_, Lx_);
145 ldl_permt(dimMat, reinterpret_cast<double*>(&x[0]), Y_, P_);
146 // this is a direct solver
147 res.iterations = 1;
148 res.converged = true;
149 }
150
152 virtual void apply(domain_type& x, range_type& b, double reduction, InverseOperatorResult& res)
153 {
154 DUNE_UNUSED_PARAMETER(reduction);
155 apply(x,b,res);
156 }
157
163 void apply(T* x, T* b)
164 {
165 const int dimMat(ldlMatrix_.N());
166 ldl_perm(dimMat, Y_, b, P_);
167 ldl_lsolve(dimMat, Y_, Lp_, Li_, Lx_);
168 ldl_dsolve(dimMat, Y_, D_);
169 ldl_ltsolve(dimMat, Y_, Lp_, Li_, Lx_);
170 ldl_permt(dimMat, x, Y_, P_);
171 }
172
173 void setOption(unsigned int option, double value)
174 {
175 DUNE_UNUSED_PARAMETER(option);
177 }
178
180 void setMatrix(const Matrix& matrix)
181 {
182 if ((ldlMatrix_.N() + ldlMatrix_.M() > 0) || matrixIsLoaded_)
183 free();
184 ldlMatrix_ = matrix;
185 decompose();
186 }
187
188 template<class S>
189 void setSubMatrix(const Matrix& matrix, const S& rowIndexSet)
190 {
191 if ((ldlMatrix_.N() + ldlMatrix_.M() > 0) || matrixIsLoaded_)
192 free();
193 ldlMatrix_.setMatrix(matrix,rowIndexSet);
194 decompose();
195 }
196
201 inline void setVerbosity(int v)
202 {
203 verbose_=v;
204 }
205
211 {
212 return ldlMatrix_;
213 }
214
219 void free()
220 {
221 delete [] D_;
222 delete [] Y_;
223 delete [] Lp_;
224 delete [] Lx_;
225 delete [] Li_;
226 delete [] P_;
227 delete [] Pinv_;
228 ldlMatrix_.free();
229 matrixIsLoaded_ = false;
230 }
231
233 inline const char* name()
234 {
235 return "LDL";
236 }
237
242 inline double* getD()
243 {
244 return D_;
245 }
246
251 inline int* getLp()
252 {
253 return Lp_;
254 }
255
260 inline int* getLi()
261 {
262 return Li_;
263 }
264
269 inline double* getLx()
270 {
271 return Lx_;
272 }
273
274 private:
275 template<class M,class X, class TM, class TD, class T1>
276 friend class SeqOverlappingSchwarz;
277
278 friend struct SeqOverlappingSchwarzAssemblerHelper<LDL<Matrix>,true>;
279
281 void decompose()
282 {
283 // allocate vectors
284 const int dimMat(ldlMatrix_.N());
285 D_ = new double [dimMat];
286 Y_ = new double [dimMat];
287 Lp_ = new int [dimMat + 1];
288 Parent_ = new int [dimMat];
289 Lnz_ = new int [dimMat];
290 Flag_ = new int [dimMat];
291 Pattern_ = new int [dimMat];
292 P_ = new int [dimMat];
293 Pinv_ = new int [dimMat];
294
295 double Info [AMD_INFO];
296 if(amd_order (dimMat, ldlMatrix_.getColStart(), ldlMatrix_.getRowIndex(), P_, (double *) NULL, Info) < AMD_OK)
297 DUNE_THROW(InvalidStateException,"Error: AMD failed!");
298 if(verbose_ > 0)
299 amd_info (Info);
300 // compute the symbolic factorisation
301 ldl_symbolic(dimMat, ldlMatrix_.getColStart(), ldlMatrix_.getRowIndex(), Lp_, Parent_, Lnz_, Flag_, P_, Pinv_);
302 // initialise those entries of additionalVectors_ whose dimension is known only now
303 Lx_ = new double [Lp_[dimMat]];
304 Li_ = new int [Lp_[dimMat]];
305 // compute the numeric factorisation
306 const int rank(ldl_numeric(dimMat, ldlMatrix_.getColStart(), ldlMatrix_.getRowIndex(), ldlMatrix_.getValues(),
307 Lp_, Parent_, Lnz_, Li_, Lx_, D_, Y_, Pattern_, Flag_, P_, Pinv_));
308 // free temporary vectors
309 delete [] Flag_;
310 delete [] Pattern_;
311 delete [] Parent_;
312 delete [] Lnz_;
313
314 if(rank!=dimMat)
315 DUNE_THROW(InvalidStateException,"Error: LDL factorisation failed!");
316 }
317
318 LDLMatrix ldlMatrix_;
319 bool matrixIsLoaded_;
320 int verbose_;
321 int* Lp_;
322 int* Parent_;
323 int* Lnz_;
324 int* Flag_;
325 int* Pattern_;
326 int* P_;
327 int* Pinv_;
328 double* D_;
329 double* Y_;
330 double* Lx_;
331 int* Li_;
332 };
333
334 template<typename T, typename A, int n, int m>
335 struct IsDirectSolver<LDL<BCRSMatrix<FieldMatrix<T,n,m>,A> > >
336 {
337 enum {value = true};
338 };
339
340 template<typename T, typename A, int n, int m>
341 struct StoresColumnCompressed<LDL<BCRSMatrix<FieldMatrix<T,n,m>,A> > >
342 {
343 enum {value = true};
344 };
345
346 struct LDLCreator {
347 template<class F> struct isValidBlock : std::false_type{};
348 template<int k> struct isValidBlock<FieldVector<double,k>> : std::true_type{};
349
350 template<typename TL, typename M>
351 std::shared_ptr<Dune::InverseOperator<typename Dune::TypeListElement<1, TL>::type,
352 typename Dune::TypeListElement<2, TL>::type>>
353 operator() (TL /*tl*/, const M& mat, const Dune::ParameterTree& config,
354 std::enable_if_t<
355 isValidBlock<typename Dune::TypeListElement<1, TL>::type::block_type>::value,int> = 0) const
356 {
357 int verbose = config.get("verbose", 0);
358 return std::make_shared<Dune::LDL<M>>(mat,verbose);
359 }
360
361 // second version with SFINAE to validate the template parameters of LDL
362 template<typename TL, typename M>
363 std::shared_ptr<Dune::InverseOperator<typename Dune::TypeListElement<1, TL>::type,
364 typename Dune::TypeListElement<2, TL>::type>>
365 operator() (TL /*tl*/, const M& /*mat*/, const Dune::ParameterTree& /*config*/,
366 std::enable_if_t<
367 !isValidBlock<typename Dune::TypeListElement<1, TL>::type::block_type>::value,int> = 0) const
368 {
369 DUNE_THROW(UnsupportedType,
370 "Unsupported Type in LDL (only double and std::complex<double> supported)");
371 }
372 };
373 DUNE_REGISTER_DIRECT_SOLVER("ldl", Dune::LDLCreator());
374
375} // end namespace Dune
376
377
378#endif //HAVE_SUITESPARSE_LDL
379#endif //DUNE_ISTL_LDL_HH
A sparse block matrix with compressed row storage.
Definition: bcrsmatrix.hh:426
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
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:97
Use the LDL package to directly solve linear systems – empty default class.
Definition: ldl.hh:55
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
Sequential overlapping Schwarz preconditioner.
Definition: overlappingschwarz.hh:752
A few common exception classes.
Implementations of the inverse operator interface.
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
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: ldl.hh:84
int * getLp()
Get factorization Lp.
Definition: ldl.hh:251
LDLMatrix & getInternalMatrix()
Return the column compress matrix.
Definition: ldl.hh:210
double * getLx()
Get factorization Lx.
Definition: ldl.hh:269
void setVerbosity(int v)
Sets the verbosity level for the solver.
Definition: ldl.hh:201
void apply(T *x, T *b)
Additional apply method with c-arrays in analogy to superlu.
Definition: ldl.hh:163
virtual void apply(domain_type &x, range_type &b, double reduction, InverseOperatorResult &res)
apply inverse operator, with given convergence criteria.
Definition: ldl.hh:152
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: ldl.hh:86
virtual ~LDL()
Default constructor.
Definition: ldl.hh:131
void free()
Free allocated space.
Definition: ldl.hh:219
virtual void apply(domain_type &x, range_type &b, InverseOperatorResult &res)
Apply inverse operator,.
Definition: ldl.hh:138
int * getLi()
Get factorization Li.
Definition: ldl.hh:260
void setMatrix(const Matrix &matrix)
Initialize data from given matrix.
Definition: ldl.hh:180
LDL(const Matrix &matrix, int verbose, bool)
Constructor for compatibility with SuperLU standard constructor.
Definition: ldl.hh:119
double * getD()
Get factorization diagonal matrix D.
Definition: ldl.hh:242
Dune::ColCompMatrix< Matrix > LDLMatrix
The corresponding SuperLU Matrix type.
Definition: ldl.hh:80
virtual SolverCategory::Category category() const
Category of the solver (see SolverCategory::Category)
Definition: ldl.hh:89
const char * name()
Get method name.
Definition: ldl.hh:233
LDL()
Default constructor.
Definition: ldl.hh:127
LDL(const Matrix &matrix, int verbose=0)
Construct a solver object from a BCRSMatrix.
Definition: ldl.hh:103
Dune namespace.
Definition: alignedallocator.hh:14
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)