3 #ifndef DUNE_ISTL_LDL_HH
4 #define DUNE_ISTL_LDL_HH
6 #if HAVE_SUITESPARSE_LDL || defined DOXYGEN
10 #include <type_traits>
20 #include <dune/common/exceptions.hh>
21 #include <dune/common/unused.hh>
41 template<
class M,
class T,
class TM,
class TD,
class TA>
42 class SeqOverlappingSchwarz;
44 template<
class T,
bool tag>
45 struct SeqOverlappingSchwarzAssemblerHelper;
53 template<
class Matrix>
70 template<
typename T,
typename A,
int n,
int m>
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> > > >
91 return SolverCategory::Category::sequential;
103 LDL(
const Matrix& matrix,
int verbose=0) : matrixIsLoaded_(false), verbose_(verbose)
106 static_assert(std::is_same<T,double>::value,
"Unsupported Type in LDL (only double supported)");
119 LDL(
const Matrix& matrix,
int verbose,
bool) : matrixIsLoaded_(false), verbose_(verbose)
122 static_assert(std::is_same<T,double>::value,
"Unsupported Type in LDL (only double supported)");
127 LDL() : matrixIsLoaded_(false), verbose_(0)
133 if ((ldlMatrix_.N() + ldlMatrix_.M() > 0) || matrixIsLoaded_)
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_);
154 DUNE_UNUSED_PARAMETER(reduction);
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_);
175 DUNE_UNUSED_PARAMETER(option);
176 DUNE_UNUSED_PARAMETER(value);
182 if ((ldlMatrix_.N() + ldlMatrix_.M() > 0) || matrixIsLoaded_)
191 if ((ldlMatrix_.N() + ldlMatrix_.M() > 0) || matrixIsLoaded_)
193 ldlMatrix_.setMatrix(matrix,rowIndexSet);
229 matrixIsLoaded_ =
false;
275 template<
class M,
class X,
class TM,
class TD,
class T1>
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];
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!");
301 ldl_symbolic(dimMat, ldlMatrix_.getColStart(), ldlMatrix_.getRowIndex(), Lp_, Parent_, Lnz_, Flag_, P_, Pinv_);
303 Lx_ =
new double [Lp_[dimMat]];
304 Li_ =
new int [Lp_[dimMat]];
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_));
315 DUNE_THROW(InvalidStateException,
"Error: LDL factorisation failed!");
318 LDLMatrix ldlMatrix_;
319 bool matrixIsLoaded_;
334 template<
typename T,
typename A,
int n,
int m>
340 template<
typename T,
typename A,
int n,
int m>
348 template<
int k>
struct isValidBlock<FieldVector<double,k>> : std::true_type{};
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>>
355 isValidBlock<
typename Dune::TypeListElement<1, TL>::type::block_type>::value,
int> = 0)
const
357 int verbose = config.get(
"verbose", 0);
358 return std::make_shared<Dune::LDL<M>>(
mat,verbose);
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>>
367 !
isValidBlock<
typename Dune::TypeListElement<1, TL>::type::block_type>::value,
int> = 0)
const
370 "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:84
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:353
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
double * getLx()
Get factorization Lx.
Definition: ldl.hh:269
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
void setMatrix(const Matrix &matrix)
Initialize data from given matrix.
Definition: ldl.hh:180
const char * name()
Get method name.
Definition: ldl.hh:233
void setOption(unsigned int option, double value)
Definition: ldl.hh:173
int * getLi()
Get factorization Li.
Definition: ldl.hh:260
LDL(const Matrix &matrix, int verbose, bool)
Constructor for compatibility with SuperLU standard constructor.
Definition: ldl.hh:119
Dune::ColCompMatrix< Matrix > LDLMatrix
The corresponding SuperLU Matrix type.
Definition: ldl.hh:80
int * getLp()
Get factorization Lp.
Definition: ldl.hh:251
virtual SolverCategory::Category category() const
Category of the solver (see SolverCategory::Category)
Definition: ldl.hh:89
void setSubMatrix(const Matrix &matrix, const S &rowIndexSet)
Definition: ldl.hh:189
LDL()
Default constructor.
Definition: ldl.hh:127
double * getD()
Get factorization diagonal matrix D.
Definition: ldl.hh:242
LDLMatrix & getInternalMatrix()
Return the column compress matrix.
Definition: ldl.hh:210
LDL(const Matrix &matrix, int verbose=0)
Construct a solver object from a BCRSMatrix.
Definition: ldl.hh:103
Matrix & mat
Definition: matrixmatrix.hh:345
Definition: allocator.hh:7
DUNE_REGISTER_DIRECT_SOLVER("cholmod", Dune::CholmodCreator())
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
Sequential overlapping Schwarz preconditioner.
Definition: overlappingschwarz.hh:752
Definition: overlappingschwarz.hh:691
Use the LDL package to directly solve linear systems – empty default class.
Definition: ldl.hh:55
Definition: matrixutils.hh:26
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
Definition: solverfactory.hh:124
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