3#ifndef DUNE_ISTL_LDL_HH
4#define DUNE_ISTL_LDL_HH
6#if HAVE_SUITESPARSE_LDL || defined DOXYGEN
19#include <dune/common/exceptions.hh>
20#include <dune/common/unused.hh>
39 template<
class M,
class T,
class TM,
class TD,
class TA>
40 class SeqOverlappingSchwarz;
42 template<
class T,
bool tag>
43 struct SeqOverlappingSchwarzAssemblerHelper;
51 template<
class Matrix>
68 template<
typename T,
typename A,
int n,
int m>
70 :
public InverseOperator<BlockVector<FieldVector<T,m>, typename A::template rebind<FieldVector<T,m> >::other>,
71 BlockVector<FieldVector<T,n>, typename A::template rebind<FieldVector<T,n> >::other> >
95 LDL(
const Matrix& matrix,
int verbose=0) : matrixIsLoaded_(
false), verbose_(verbose)
98 static_assert(std::is_same<T,double>::value,
"Unsupported Type in LDL (only double supported)");
111 LDL(
const Matrix& matrix,
int verbose,
bool) : matrixIsLoaded_(
false), verbose_(verbose)
114 static_assert(std::is_same<T,double>::value,
"Unsupported Type in LDL (only double supported)");
125 if ((ldlMatrix_.N() + ldlMatrix_.M() > 0) || matrixIsLoaded_)
132 const int dimMat(ldlMatrix_.N());
140 res.converged =
true;
157 const int dimMat(ldlMatrix_.N());
174 if ((ldlMatrix_.N() + ldlMatrix_.M() > 0) || matrixIsLoaded_)
183 if ((ldlMatrix_.N() + ldlMatrix_.M() > 0) || matrixIsLoaded_)
185 ldlMatrix_.setMatrix(matrix,rowIndexSet);
221 matrixIsLoaded_ =
false;
267 template<
class M,
class X,
class TM,
class TD,
class T1>
276 const int dimMat(ldlMatrix_.N());
279 Lp_ =
new int [
dimMat + 1];
280 Parent_ =
new int [
dimMat];
283 Pattern_ =
new int [
dimMat];
293 ldl_symbolic(
dimMat, ldlMatrix_.getColStart(), ldlMatrix_.getRowIndex(), Lp_, Parent_, Lnz_, Flag_, P_, Pinv_);
295 Lx_ =
new double [Lp_[
dimMat]];
296 Li_ =
new int [Lp_[
dimMat]];
298 const int rank(
ldl_numeric(
dimMat, ldlMatrix_.getColStart(), ldlMatrix_.getRowIndex(), ldlMatrix_.getValues(),
299 Lp_, Parent_, Lnz_, Li_, Lx_, D_, Y_, Pattern_, Flag_, P_, Pinv_));
310 LDLMatrix ldlMatrix_;
311 bool matrixIsLoaded_;
326 template<
typename T,
typename A,
int n,
int m>
332 template<
typename T,
typename A,
int n,
int m>
Implementations of the inverse operator interface.
Templates characterizing the type of a solver.
Dune::BCRSMatrix< FieldMatrix< T, n, m >, A > Matrix
The matrix type.
Definition ldl.hh:75
int * getLp()
Get factorization Lp.
Definition ldl.hh:243
LDLMatrix & getInternalMatrix()
Return the column compress matrix.
Definition ldl.hh:202
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:84
double * getLx()
Get factorization Lx.
Definition ldl.hh:261
void setVerbosity(int v)
Sets the verbosity level for the solver.
Definition ldl.hh:193
void apply(T *x, T *b)
Additional apply method with c-arrays in analogy to superlu.
Definition ldl.hh:155
ColCompMatrixInitializer< BCRSMatrix< FieldMatrix< T, n, m >, A > > MatrixInitializer
Type of an associated initializer class.
Definition ldl.hh:80
virtual void apply(domain_type &x, range_type &b, double reduction, InverseOperatorResult &res)
apply inverse operator, with given convergence criteria.
Definition ldl.hh:144
virtual ~LDL()
Default constructor.
Definition ldl.hh:123
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:82
void free()
Free allocated space.
Definition ldl.hh:211
virtual void apply(domain_type &x, range_type &b, InverseOperatorResult &res)
Apply inverse operator,.
Definition ldl.hh:130
int * getLi()
Get factorization Li.
Definition ldl.hh:252
Dune::BCRSMatrix< FieldMatrix< T, n, m >, A > matrix_type
Definition ldl.hh:76
void setMatrix(const Matrix &matrix)
Initialize data from given matrix.
Definition ldl.hh:172
void setOption(unsigned int option, double value)
Definition ldl.hh:165
LDL(const Matrix &matrix, int verbose, bool)
Constructor for compatibility with SuperLU standard constructor.
Definition ldl.hh:111
double * getD()
Get factorization diagonal matrix D.
Definition ldl.hh:234
Dune::ColCompMatrix< Matrix > LDLMatrix
The corresponding SuperLU Matrix type.
Definition ldl.hh:78
const char * name()
Get method name.
Definition ldl.hh:225
void setSubMatrix(const Matrix &matrix, const S &rowIndexSet)
Definition ldl.hh:181
LDL()
Default constructor.
Definition ldl.hh:119
LDL(const Matrix &matrix, int verbose=0)
Construct a solver object from a BCRSMatrix.
Definition ldl.hh:95
Definition basearray.hh:19
Statistics about compression achieved in implicit mode.
Definition bcrsmatrix.hh:81
A sparse block matrix with compressed row storage.
Definition bcrsmatrix.hh:412
A vector of blocks with memory management.
Definition bvector.hh:309
Sequential overlapping Schwarz preconditioner.
Definition overlappingschwarz.hh:742
Definition overlappingschwarz.hh:683
Use the LDL package to directly solve linear systems – empty default class.
Definition ldl.hh:53
Definition matrixutils.hh:25
Statistics about the application of an inverse operator.
Definition solver.hh:32
Abstract base class for all solvers.
Definition solver.hh:79
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