DUNE PDELab (git)
spqr.hh
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright © DUNE Project contributors, see file LICENSE.md in module root
65 : public InverseOperator<BlockVector<FieldVector<T,m>, typename std::allocator_traits<A>::template rebind_alloc<FieldVector<T,m> > >,
66 BlockVector<FieldVector<T,n>, typename std::allocator_traits<A>::template rebind_alloc<FieldVector<T,n> > > >
75 typedef ISTL::Impl::BCCSMatrixInitializer<BCRSMatrix<FieldMatrix<T,n,m>,A>, int> MatrixInitializer;
77 typedef Dune::BlockVector<FieldVector<T,m>, typename std::allocator_traits<A>::template rebind_alloc<FieldVector<T,m> > > domain_type;
79 typedef Dune::BlockVector<FieldVector<T,n>, typename std::allocator_traits<A>::template rebind_alloc<FieldVector<T,n> > > range_type;
116 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
140 static_assert((std::is_same<T,double>::value) || (std::is_same<T,std::complex<double> >::value),
191 void apply (domain_type& x, range_type& b, [[maybe_unused]] double reduction, InverseOperatorResult& res) override
229 copyToBCCSMatrix(initializer, ISTL::Impl::MatrixRowSubset<Matrix,std::set<std::size_t> >(matrix,rowIndexSet));
Default exception if a function was called while the object is not in a valid state for that function...
Definition: exceptions.hh:281
std::string get(const std::string &key, const std::string &defaultValue) const
get value as string
Definition: parametertree.cc:188
Use the SPQR package to directly solve linear systems – empty default class.
Definition: spqr.hh:48
Sequential overlapping Schwarz preconditioner.
Definition: overlappingschwarz.hh:755
A few common exception classes.
Implementations of the inverse operator interface.
virtual ~SPQR()
Destructor.
Definition: spqr.hh:147
SPQR(const Matrix &matrix, int verbose, bool)
Constructor for compatibility with SuperLU standard constructor.
Definition: spqr.hh:113
SPQRMatrix & getInternalMatrix()
Return the column coppressed matrix.
Definition: spqr.hh:256
void setMatrix(const Matrix &matrix)
Initialize data from given matrix.
Definition: spqr.hh:200
void apply(domain_type &x, range_type &b, InverseOperatorResult &res) override
Apply inverse operator,.
Definition: spqr.hh:155
ISTL::Impl::BCCSMatrixInitializer< BCRSMatrix< FieldMatrix< T, n, m >, A >, int > MatrixInitializer
Type of an associated initializer class.
Definition: spqr.hh:75
SolverCategory::Category category() const override
Category of the solver (see SolverCategory::Category)
Definition: spqr.hh:82
SPQR()
Default constructor.
Definition: spqr.hh:137
const char * name()
Get method name.
Definition: spqr.hh:275
SuiteSparseQR_factorization< T > * getFactorization()
Return the matrix factorization.
Definition: spqr.hh:247
void setVerbosity(int v)
Sets the verbosity level for the solver.
Definition: spqr.hh:238
ISTL::Impl::BCCSMatrix< T, int > SPQRMatrix
The corresponding SuperLU Matrix type.
Definition: spqr.hh:73
void free()
Free allocated space.
Definition: spqr.hh:265
SPQR(const Matrix &matrix, const ParameterTree &config)
Constructs the SPQR solver.
Definition: spqr.hh:132
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:79
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:77
SPQR(const Matrix &matrix, int verbose=0)
Construct a solver object from a BCRSMatrix.
Definition: spqr.hh:95
void apply(domain_type &x, range_type &b, double reduction, InverseOperatorResult &res) override
apply inverse operator, with given convergence criteria.
Definition: spqr.hh:191
constexpr auto get(std::integer_sequence< T, II... >, std::integral_constant< std::size_t, pos >={})
Return the entry at position pos of the given sequence.
Definition: integersequence.hh:22
Templates characterizing the type of a solver.
Statistics about the application of an inverse operator.
Definition: solver.hh:50
bool converged
True if convergence criterion has been met.
Definition: solver.hh:75
|
Legal Statements / Impressum |
Hosted by TU Dresden |
generated with Hugo v0.111.3
(Nov 12, 23:30, 2024)