KASKADE 7 development version
directPreconditioner.hh
Go to the documentation of this file.
1/* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
2/* */
3/* This file is part of the library KASKADE 7 */
4/* https://www.zib.de/research/projects/kaskade7-finite-element-toolbox */
5/* */
6/* Copyright (C) 2002-2016 Zuse Institute Berlin */
7/* */
8/* KASKADE 7 is distributed under the terms of the ZIB Academic License. */
9/* see $KASKADE/academic.txt */
10/* */
11/* * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * */
12
13#ifndef DIRECT_PRECONDITIONER_HH
14#define DIRECT_PRECONDITIONER_HH
15
16#include <memory>
17
18#include <dune/istl/preconditioners.hh>
19#include <dune/istl/solvercategory.hh>
20
21#include "linalg/crsutil.hh"
25
26namespace Kaskade
27{
29 // forward declarations
30 template <class,class> class Factorization;
31 template <class,class> class NumaBCRSMatrix;
33
43 template <class Op>
44 class DirectPreconditioner: public Dune::Preconditioner<typename Op::Domain, typename Op::Range>
45 {
46 typedef typename Op::Domain Domain;
47 typedef typename Op::Domain Range;
48 typedef typename GetScalar<Domain>::type Scalar;
49
50 public:
51 static int const category = Dune::SolverCategory::sequential;
52
60 : factorization(getFactorization(directType,op.template get<MatrixAsTriplet<Scalar> >(),typename Factorization<Scalar>::Options(property)))
61 {}
62
63 void pre(Domain&, Range&) {}
64 void post (Domain&) {}
65
66 void apply (Domain& x, Range const& y)
67 {
68 std::vector<Scalar> b(y.dim());
69 vectorToSequence(y,begin(b));
70 factorization->solve(b);
71 vectorFromSequence(x,begin(b));
72 }
73
74 private:
75 std::unique_ptr<Factorization<Scalar>> factorization;
76 };
77
78
79
80} // namespace Kaskade
81#endif
void apply(Domain &x, Range const &y)
DirectPreconditioner(Op const &op, DirectType directType=DirectType::UMFPACK, MatrixProperties property=MatrixProperties::SYMMETRIC)
Abstract base class for matrix factorizations.
std::unique_ptr< Factorization< Scalar > > getFactorization(DirectType directType, MatrixAsTriplet< Scalar, Index > const &A, FactorizationOptions options)
Creates a factorization of the given triplet matrix.
DirectType
Available direct solvers for linear equation systems.
Definition: enums.hh:35
@ UMFPACK
UMFPack from SuiteSparse, using 32 bit integer indices.
MatrixProperties
Characterizations of sparse matrix properties.
Definition: enums.hh:76
OutIter vectorToSequence(double x, OutIter iter)
Definition: crsutil.hh:30
InIter vectorFromSequence(FunctionSpaceElement< Space, m > &v, InIter in)
!std ::is_same< decltype(hasFieldType< Type >(0)), TypeNotFound >::value ::type type