3#ifndef DUNE_ISTL_LDL_HH
4#define DUNE_ISTL_LDL_HH
6#if HAVE_SUITESPARSE_LDL || defined DOXYGEN
20#include <dune/common/exceptions.hh>
40 template<
class M,
class T,
class TM,
class TD,
class TA>
41 class SeqOverlappingSchwarz;
43 template<
class T,
bool tag>
44 struct SeqOverlappingSchwarzAssemblerHelper;
52 template<
class Matrix>
69 template<
typename T,
typename A,
int n,
int m>
71 :
public InverseOperator<BlockVector<FieldVector<T,m>, typename std::allocator_traits<A>::template rebind_alloc<FieldVector<T,m> > >,
72 BlockVector<FieldVector<T,n>, typename std::allocator_traits<A>::template rebind_alloc<FieldVector<T,n> > > >
79 typedef Dune::ISTL::Impl::BCCSMatrix<T,int>
LDLMatrix;
81 typedef ISTL::Impl::BCCSMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A>,
int>
MatrixInitializer;
102 LDL(
const Matrix& matrix,
int verbose=0) : matrixIsLoaded_(false), verbose_(verbose)
105 static_assert(std::is_same<T,double>::value,
"Unsupported Type in LDL (only double supported)");
118 LDL(
const Matrix& matrix,
int verbose,
bool) : matrixIsLoaded_(false), verbose_(verbose)
121 static_assert(std::is_same<T,double>::value,
"Unsupported Type in LDL (only double supported)");
135 :
LDL(matrix, config.
get<int>(
"verbose", 0))
139 LDL() : matrixIsLoaded_(false), verbose_(0)
145 if ((ldlMatrix_.N() + ldlMatrix_.M() > 0) || matrixIsLoaded_)
152 const int dimMat(ldlMatrix_.N());
153 ldl_perm(dimMat, Y_,
reinterpret_cast<double*
>(&b[0]), P_);
154 ldl_lsolve(dimMat, Y_, Lp_, Li_, Lx_);
155 ldl_dsolve(dimMat, Y_, D_);
156 ldl_ltsolve(dimMat, Y_, Lp_, Li_, Lx_);
157 ldl_permt(dimMat,
reinterpret_cast<double*
>(&x[0]), Y_, P_);
176 const int dimMat(ldlMatrix_.N());
177 ldl_perm(dimMat, Y_, b, P_);
178 ldl_lsolve(dimMat, Y_, Lp_, Li_, Lx_);
179 ldl_dsolve(dimMat, Y_, D_);
180 ldl_ltsolve(dimMat, Y_, Lp_, Li_, Lx_);
181 ldl_permt(dimMat, x, Y_, P_);
184 void setOption([[maybe_unused]]
unsigned int option, [[maybe_unused]]
double value)
190 if ((ldlMatrix_.N() + ldlMatrix_.M() > 0) || matrixIsLoaded_)
193 if (ldlMatrix_.N() + ldlMatrix_.M() + ldlMatrix_.nonzeroes() != 0)
197 ISTL::Impl::BCCSMatrixInitializer<Matrix, int> initializer(ldlMatrix_);
199 copyToBCCSMatrix(initializer, matrix);
207 if ((ldlMatrix_.N() + ldlMatrix_.M() > 0) || matrixIsLoaded_)
210 if (ldlMatrix_.N() + ldlMatrix_.M() + ldlMatrix_.nonzeroes() != 0)
215 ISTL::Impl::BCCSMatrixInitializer<Matrix, int> initializer(ldlMatrix_);
217 copyToBCCSMatrix(initializer, ISTL::Impl::MatrixRowSubset<
Matrix,std::set<std::size_t> >(matrix,rowIndexSet));
254 matrixIsLoaded_ =
false;
300 template<
class M,
class X,
class TM,
class TD,
class T1>
309 const int dimMat(ldlMatrix_.N());
310 D_ =
new double [dimMat];
311 Y_ =
new double [dimMat];
312 Lp_ =
new int [dimMat + 1];
313 Parent_ =
new int [dimMat];
314 Lnz_ =
new int [dimMat];
315 Flag_ =
new int [dimMat];
316 Pattern_ =
new int [dimMat];
317 P_ =
new int [dimMat];
318 Pinv_ =
new int [dimMat];
320 double Info [AMD_INFO];
321 if(amd_order (dimMat, ldlMatrix_.getColStart(), ldlMatrix_.getRowIndex(), P_, (
double *) NULL, Info) < AMD_OK)
322 DUNE_THROW(InvalidStateException,
"Error: AMD failed!");
326 ldl_symbolic(dimMat, ldlMatrix_.getColStart(), ldlMatrix_.getRowIndex(), Lp_, Parent_, Lnz_, Flag_, P_, Pinv_);
328 Lx_ =
new double [Lp_[dimMat]];
329 Li_ =
new int [Lp_[dimMat]];
331 const int rank(ldl_numeric(dimMat, ldlMatrix_.getColStart(), ldlMatrix_.getRowIndex(), ldlMatrix_.getValues(),
332 Lp_, Parent_, Lnz_, Li_, Lx_, D_, Y_, Pattern_, Flag_, P_, Pinv_));
340 DUNE_THROW(InvalidStateException,
"Error: LDL factorisation failed!");
343 LDLMatrix ldlMatrix_;
344 bool matrixIsLoaded_;
359 template<
typename T,
typename A,
int n,
int m>
365 template<
typename T,
typename A,
int n,
int m>
373 template<
int k>
struct isValidBlock<FieldVector<double,k>> : std::true_type{};
375 template<
typename TL,
typename M>
376 std::shared_ptr<Dune::InverseOperator<typename Dune::TypeListElement<1, TL>::type,
377 typename Dune::TypeListElement<2, TL>::type>>
380 isValidBlock<
typename Dune::TypeListElement<1, TL>::type::block_type>::value,
int> = 0)
const
382 int verbose = config.get(
"verbose", 0);
383 return std::make_shared<Dune::LDL<M>>(
mat,verbose);
387 template<
typename TL,
typename M>
388 std::shared_ptr<Dune::InverseOperator<typename Dune::TypeListElement<1, TL>::type,
389 typename Dune::TypeListElement<2, TL>::type>>
392 !
isValidBlock<
typename Dune::TypeListElement<1, TL>::type::block_type>::value,
int> = 0)
const
395 "Unsupported Type in LDL (only double and std::complex<double> supported)");
Templates characterizing the type of a solver.
Implementations of the inverse operator interface.
#define DUNE_REGISTER_DIRECT_SOLVER(name,...)
Definition: solverregistry.hh:11
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:83
Dune::BCRSMatrix< FieldMatrix< T, n, m >, A > Matrix
The matrix type.
Definition: ldl.hh:76
int * getLp()
Get factorization Lp.
Definition: ldl.hh:276
LDLMatrix & getInternalMatrix()
Return the column compress matrix.
Definition: ldl.hh:235
Dune::ISTL::Impl::BCCSMatrix< T, int > LDLMatrix
The corresponding SuperLU Matrix type.
Definition: ldl.hh:79
double * getLx()
Get factorization Lx.
Definition: ldl.hh:294
void setVerbosity(int v)
Sets the verbosity level for the solver.
Definition: ldl.hh:226
void apply(T *x, T *b)
Additional apply method with c-arrays in analogy to superlu.
Definition: ldl.hh:174
virtual void apply(domain_type &x, range_type &b, double reduction, InverseOperatorResult &res)
apply inverse operator, with given convergence criteria.
Definition: ldl.hh:164
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:85
virtual ~LDL()
Default constructor.
Definition: ldl.hh:143
void free()
Free allocated space.
Definition: ldl.hh:244
virtual void apply(domain_type &x, range_type &b, InverseOperatorResult &res)
Apply inverse operator,.
Definition: ldl.hh:150
int * getLi()
Get factorization Li.
Definition: ldl.hh:285
Dune::BCRSMatrix< FieldMatrix< T, n, m >, A > matrix_type
Definition: ldl.hh:77
void setMatrix(const Matrix &matrix)
Initialize data from given matrix.
Definition: ldl.hh:188
void setOption(unsigned int option, double value)
Definition: ldl.hh:184
LDL(const Matrix &matrix, int verbose, bool)
Constructor for compatibility with SuperLU standard constructor.
Definition: ldl.hh:118
double * getD()
Get factorization diagonal matrix D.
Definition: ldl.hh:267
LDL(const Matrix &matrix, const ParameterTree &config)
Constructs the LDL solver.
Definition: ldl.hh:134
virtual SolverCategory::Category category() const
Category of the solver (see SolverCategory::Category)
Definition: ldl.hh:88
ISTL::Impl::BCCSMatrixInitializer< BCRSMatrix< FieldMatrix< T, n, m >, A >, int > MatrixInitializer
Type of an associated initializer class.
Definition: ldl.hh:81
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:378
const char * name()
Get method name.
Definition: ldl.hh:258
void setSubMatrix(const Matrix &matrix, const S &rowIndexSet)
Definition: ldl.hh:205
LDL()
Default constructor.
Definition: ldl.hh:139
LDL(const Matrix &matrix, int verbose=0)
Construct a solver object from a BCRSMatrix.
Definition: ldl.hh:102
Matrix & mat
Definition: matrixmatrix.hh:345
Definition: allocator.hh:9
PropertyMapTypeSelector< Amg::VertexVisitedTag, Amg::PropertiesGraph< G, Amg::VertexProperties, EP, VM, EM > >::Type get(const Amg::VertexVisitedTag &tag, Amg::PropertiesGraph< G, Amg::VertexProperties, EP, VM, EM > &graph)
Definition: dependency.hh:291
Definition: matrixutils.hh:209
A sparse block matrix with compressed row storage.
Definition: bcrsmatrix.hh:464
size_type M() const
number of columns (counted in blocks)
Definition: bcrsmatrix.hh:1976
size_type N() const
number of rows (counted in blocks)
Definition: bcrsmatrix.hh:1970
A vector of blocks with memory management.
Definition: bvector.hh:393
Sequential overlapping Schwarz preconditioner.
Definition: overlappingschwarz.hh:753
Definition: overlappingschwarz.hh:692
Use the LDL package to directly solve linear systems – empty default class.
Definition: ldl.hh:54
Definition: matrixutils.hh:25
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
Abstract base class for all solvers.
Definition: solver.hh:97
Category
Definition: solvercategory.hh:21
@ sequential
Category for sequential solvers.
Definition: solvercategory.hh:23
Definition: solverregistry.hh:75
Definition: solvertype.hh:14
@ value
Whether this is a direct solver.
Definition: solvertype.hh:22
Definition: solvertype.hh:28
@ value
whether the solver internally uses column compressed storage
Definition: solvertype.hh:34