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