Dune Core Modules (2.7.1)
•All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
spqr.hh
Go to the documentation of this file.
64 : public InverseOperator<BlockVector<FieldVector<T,m>, typename std::allocator_traits<A>::template rebind_alloc<FieldVector<T,m> > >,
65 BlockVector<FieldVector<T,n>, typename std::allocator_traits<A>::template rebind_alloc<FieldVector<T,n> > > >
76 typedef Dune::BlockVector<FieldVector<T,m>, typename std::allocator_traits<A>::template rebind_alloc<FieldVector<T,m> > > domain_type;
78 typedef Dune::BlockVector<FieldVector<T,n>, typename std::allocator_traits<A>::template rebind_alloc<FieldVector<T,n> > > range_type;
115 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
126 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
171 virtual void apply (domain_type& x, range_type& b, double reduction, InverseOperatorResult& res)
314 std::enable_if_t<!isValidBlock<typename Dune::TypeListElement<1, TL>::type::block_type>::value,int> = 0) const
322 //template<> struct SPQRCreator::isValidMatrixBlock<FieldMatrix<std::complex<double>,1,1>> : std::true_type{};
Inititializer for the ColCompMatrix as needed by OverlappingSchwarz.
Definition: colcompmatrix.hh:255
std::string get(const std::string &key, const std::string &defaultValue) const
get value as string
Definition: parametertree.cc:183
Use the SPQR package to directly solve linear systems – empty default class.
Definition: spqr.hh:47
Sequential overlapping Schwarz preconditioner.
Definition: overlappingschwarz.hh:752
A few common exception classes.
decltype(auto) apply(F &&f, ArgTuple &&args)
Apply function with arguments given as tuple.
Definition: apply.hh:46
#define DUNE_UNUSED_PARAMETER(parm)
A macro to mark intentionally unused function parameters with.
Definition: unused.hh:25
virtual ~SPQR()
Destructor.
Definition: spqr.hh:133
SPQR(const Matrix &matrix, int verbose, bool)
Constructor for compatibility with SuperLU standard constructor.
Definition: spqr.hh:112
virtual SolverCategory::Category category() const
Category of the solver (see SolverCategory::Category)
Definition: spqr.hh:81
SPQRMatrix & getInternalMatrix()
Return the column coppressed matrix.
Definition: spqr.hh:223
void setMatrix(const Matrix &matrix)
Initialize data from given matrix.
Definition: spqr.hh:184
SPQR()
Default constructor.
Definition: spqr.hh:123
const char * name()
Get method name.
Definition: spqr.hh:242
SuiteSparseQR_factorization< T > * getFactorization()
Return the matrix factorization.
Definition: spqr.hh:214
void setVerbosity(int v)
Sets the verbosity level for the solver.
Definition: spqr.hh:205
Dune::ColCompMatrix< Matrix > SPQRMatrix
The corresponding SuperLU Matrix type.
Definition: spqr.hh:72
virtual void apply(domain_type &x, range_type &b, double reduction, InverseOperatorResult &res)
apply inverse operator, with given convergence criteria.
Definition: spqr.hh:171
void free()
Free allocated space.
Definition: spqr.hh:232
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: spqr.hh:78
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: spqr.hh:76
virtual void apply(domain_type &x, range_type &b, InverseOperatorResult &res)
Apply inverse operator,.
Definition: spqr.hh:141
SPQR(const Matrix &matrix, int verbose=0)
Construct a solver object from a BCRSMatrix.
Definition: spqr.hh:94
Implementations of the inverse operator interface.
Templates characterizing the type of a solver.
Statistics about the application of an inverse operator.
Definition: solver.hh:46
bool converged
True if convergence criterion has been met.
Definition: solver.hh:71
Definition of the DUNE_UNUSED macro for the case that config.h is not available.
