5#ifndef DUNE_ISTL_LDL_HH
6#define DUNE_ISTL_LDL_HH
8#if HAVE_SUITESPARSE_LDL || defined DOXYGEN
22#include <dune/common/exceptions.hh>
42 template<
class M,
class T,
class TM,
class TD,
class TA>
43 class SeqOverlappingSchwarz;
45 template<
class T,
bool tag>
46 struct SeqOverlappingSchwarzAssemblerHelper;
54 template<
class Matrix>
71 template<
typename T,
typename A,
int n,
int m>
73 :
public InverseOperator<BlockVector<FieldVector<T,m>, typename std::allocator_traits<A>::template rebind_alloc<FieldVector<T,m> > >,
74 BlockVector<FieldVector<T,n>, typename std::allocator_traits<A>::template rebind_alloc<FieldVector<T,n> > > >
81 typedef Dune::ISTL::Impl::BCCSMatrix<T,int>
LDLMatrix;
83 typedef ISTL::Impl::BCCSMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A>,
int>
MatrixInitializer;
92 return SolverCategory::Category::sequential;
104 LDL(
const Matrix& matrix,
int verbose=0) : matrixIsLoaded_(false), verbose_(verbose)
107 static_assert(std::is_same<T,double>::value,
"Unsupported Type in LDL (only double supported)");
120 LDL(
const Matrix& matrix,
int verbose,
bool) : matrixIsLoaded_(false), verbose_(verbose)
123 static_assert(std::is_same<T,double>::value,
"Unsupported Type in LDL (only double supported)");
137 :
LDL(matrix, config.
get<int>(
"verbose", 0))
141 LDL() : matrixIsLoaded_(false), verbose_(0)
147 if ((ldlMatrix_.N() + ldlMatrix_.M() > 0) || matrixIsLoaded_)
154 const int dimMat(ldlMatrix_.N());
155 ldl_perm(dimMat, Y_,
reinterpret_cast<double*
>(&b[0]), P_);
156 ldl_lsolve(dimMat, Y_, Lp_, Li_, Lx_);
157 ldl_dsolve(dimMat, Y_, D_);
158 ldl_ltsolve(dimMat, Y_, Lp_, Li_, Lx_);
159 ldl_permt(dimMat,
reinterpret_cast<double*
>(&x[0]), Y_, P_);
178 const int dimMat(ldlMatrix_.N());
179 ldl_perm(dimMat, Y_, b, P_);
180 ldl_lsolve(dimMat, Y_, Lp_, Li_, Lx_);
181 ldl_dsolve(dimMat, Y_, D_);
182 ldl_ltsolve(dimMat, Y_, Lp_, Li_, Lx_);
183 ldl_permt(dimMat, x, Y_, P_);
186 void setOption([[maybe_unused]]
unsigned int option, [[maybe_unused]]
double value)
192 if ((ldlMatrix_.N() + ldlMatrix_.M() > 0) || matrixIsLoaded_)
195 if (ldlMatrix_.N() + ldlMatrix_.M() + ldlMatrix_.nonzeroes() != 0)
199 ISTL::Impl::BCCSMatrixInitializer<Matrix, int> initializer(ldlMatrix_);
201 copyToBCCSMatrix(initializer, matrix);
209 if ((ldlMatrix_.N() + ldlMatrix_.M() > 0) || matrixIsLoaded_)
212 if (ldlMatrix_.N() + ldlMatrix_.M() + ldlMatrix_.nonzeroes() != 0)
217 ISTL::Impl::BCCSMatrixInitializer<Matrix, int> initializer(ldlMatrix_);
219 copyToBCCSMatrix(initializer, ISTL::Impl::MatrixRowSubset<
Matrix,std::set<std::size_t> >(matrix,rowIndexSet));
256 matrixIsLoaded_ =
false;
302 template<
class M,
class X,
class TM,
class TD,
class T1>
311 const int dimMat(ldlMatrix_.N());
312 D_ =
new double [dimMat];
313 Y_ =
new double [dimMat];
314 Lp_ =
new int [dimMat + 1];
315 Parent_ =
new int [dimMat];
316 Lnz_ =
new int [dimMat];
317 Flag_ =
new int [dimMat];
318 Pattern_ =
new int [dimMat];
319 P_ =
new int [dimMat];
320 Pinv_ =
new int [dimMat];
322 double Info [AMD_INFO];
323 if(amd_order (dimMat, ldlMatrix_.getColStart(), ldlMatrix_.getRowIndex(), P_, (
double *) NULL, Info) < AMD_OK)
324 DUNE_THROW(InvalidStateException,
"Error: AMD failed!");
328 ldl_symbolic(dimMat, ldlMatrix_.getColStart(), ldlMatrix_.getRowIndex(), Lp_, Parent_, Lnz_, Flag_, P_, Pinv_);
330 Lx_ =
new double [Lp_[dimMat]];
331 Li_ =
new int [Lp_[dimMat]];
333 const int rank(ldl_numeric(dimMat, ldlMatrix_.getColStart(), ldlMatrix_.getRowIndex(), ldlMatrix_.getValues(),
334 Lp_, Parent_, Lnz_, Li_, Lx_, D_, Y_, Pattern_, Flag_, P_, Pinv_));
342 DUNE_THROW(InvalidStateException,
"Error: LDL factorisation failed!");
345 LDLMatrix ldlMatrix_;
346 bool matrixIsLoaded_;
361 template<
typename T,
typename A,
int n,
int m>
367 template<
typename T,
typename A,
int n,
int m>
375 template<
int k>
struct isValidBlock<FieldVector<double,k>> : std::true_type{};
377 template<
typename TL,
typename M>
378 std::shared_ptr<Dune::InverseOperator<typename Dune::TypeListElement<1, TL>::type,
379 typename Dune::TypeListElement<2, TL>::type>>
382 isValidBlock<
typename Dune::TypeListElement<1, TL>::type::block_type>::value,
int> = 0)
const
384 int verbose = config.get(
"verbose", 0);
385 return std::make_shared<Dune::LDL<M>>(
mat,verbose);
389 template<
typename TL,
typename M>
390 std::shared_ptr<Dune::InverseOperator<typename Dune::TypeListElement<1, TL>::type,
391 typename Dune::TypeListElement<2, TL>::type>>
394 !
isValidBlock<
typename Dune::TypeListElement<1, TL>::type::block_type>::value,
int> = 0)
const
397 "Unsupported Type in LDL (only double and std::complex<double> supported)");
Implementations of the inverse operator interface.
Templates characterizing the type of a solver.
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:85
int * getLp()
Get factorization Lp.
Definition: ldl.hh:278
LDLMatrix & getInternalMatrix()
Return the column compress matrix.
Definition: ldl.hh:237
Dune::ISTL::Impl::BCCSMatrix< T, int > LDLMatrix
The corresponding SuperLU Matrix type.
Definition: ldl.hh:81
double * getLx()
Get factorization Lx.
Definition: ldl.hh:296
DUNE_REGISTER_DIRECT_SOLVER("ldl", Dune::LDLCreator())
void setVerbosity(int v)
Sets the verbosity level for the solver.
Definition: ldl.hh:228
void apply(T *x, T *b)
Additional apply method with c-arrays in analogy to superlu.
Definition: ldl.hh:176
virtual void apply(domain_type &x, range_type &b, double reduction, InverseOperatorResult &res)
apply inverse operator, with given convergence criteria.
Definition: ldl.hh:166
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:87
virtual ~LDL()
Default constructor.
Definition: ldl.hh:145
void free()
Free allocated space.
Definition: ldl.hh:246
virtual void apply(domain_type &x, range_type &b, InverseOperatorResult &res)
Apply inverse operator,.
Definition: ldl.hh:152
int * getLi()
Get factorization Li.
Definition: ldl.hh:287
void setMatrix(const Matrix &matrix)
Initialize data from given matrix.
Definition: ldl.hh:190
void setOption(unsigned int option, double value)
Definition: ldl.hh:186
LDL(const Matrix &matrix, int verbose, bool)
Constructor for compatibility with SuperLU standard constructor.
Definition: ldl.hh:120
double * getD()
Get factorization diagonal matrix D.
Definition: ldl.hh:269
LDL(const Matrix &matrix, const ParameterTree &config)
Constructs the LDL solver.
Definition: ldl.hh:136
virtual SolverCategory::Category category() const
Category of the solver (see SolverCategory::Category)
Definition: ldl.hh:90
ISTL::Impl::BCCSMatrixInitializer< BCRSMatrix< FieldMatrix< T, n, m >, A >, int > MatrixInitializer
Type of an associated initializer class.
Definition: ldl.hh:83
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:380
const char * name()
Get method name.
Definition: ldl.hh:260
void setSubMatrix(const Matrix &matrix, const S &rowIndexSet)
Definition: ldl.hh:207
LDL()
Default constructor.
Definition: ldl.hh:141
LDL(const Matrix &matrix, int verbose=0)
Construct a solver object from a BCRSMatrix.
Definition: ldl.hh:104
Matrix & mat
Definition: matrixmatrix.hh:347
Definition: allocator.hh:11
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:293
Definition: matrixutils.hh:211
A sparse block matrix with compressed row storage.
Definition: bcrsmatrix.hh:466
size_type M() const
number of columns (counted in blocks)
Definition: bcrsmatrix.hh:1978
size_type N() const
number of rows (counted in blocks)
Definition: bcrsmatrix.hh:1972
A vector of blocks with memory management.
Definition: bvector.hh:395
Sequential overlapping Schwarz preconditioner.
Definition: overlappingschwarz.hh:755
Definition: overlappingschwarz.hh:694
Use the LDL package to directly solve linear systems – empty default class.
Definition: ldl.hh:56
Definition: matrixutils.hh:27
Statistics about the application of an inverse operator.
Definition: solver.hh:48
int iterations
Number of iterations.
Definition: solver.hh:67
bool converged
True if convergence criterion has been met.
Definition: solver.hh:73
Abstract base class for all solvers.
Definition: solver.hh:99
Category
Definition: solvercategory.hh:23
Definition: solverregistry.hh:77
Definition: solvertype.hh:16
@ value
Whether this is a direct solver.
Definition: solvertype.hh:24
Definition: solvertype.hh:30
@ value
whether the solver internally uses column compressed storage
Definition: solvertype.hh:36