Dune::Copasi
Loading...
Searching...
No Matches
inverse.hh
Go to the documentation of this file.
1#ifndef DUNE_COPASI_SOLVER_ISTL_FACTORY_INVERSE_HH
2#define DUNE_COPASI_SOLVER_ISTL_FACTORY_INVERSE_HH
3
7
8#include <dune/istl/scalarproducts.hh>
9#include <dune/istl/solver.hh>
10
11#ifdef DUNE_COPASI_DEFAULT_DIRECT_SOLVER
12#define DUNE_COPASI_DEFAULT_LINEAR_INVERSE_OPERATOR DUNE_COPASI_DEFAULT_DIRECT_SOLVER
13#else
14#define DUNE_COPASI_DEFAULT_LINEAR_INVERSE_OPERATOR DUNE_COPASI_DEFAULT_ITERATIVE_SOLVER
15#endif
16
17namespace Dune::Copasi::ISTL {
18
19namespace Impl {
20
21template<Concept::LinearOperator Op, class Alloc = std::allocator<Op>>
22std::shared_ptr<InverseOperator<typename Op::domain_type, typename Op::range_type>>
23makeInverseOperator(const std::shared_ptr<Op>& op,
24 const ParameterTree& config,
25 const Alloc& alloc = Alloc(),
26 ADLTag = {})
27{
28 assert(op);
30 if (config.get("verbosity", 1) > 0)
31 spdlog::info("Creating linear solver with type '{}'", type_name);
33 auto prec = makePreconditioner(op, config.sub("preconditioner"));
34 if constexpr (requires { op->getCommunication(); })
35 throw format_exception(NotImplemented{}, "Parallel solvers have not been implemented!");
37 using ScalarProdAlloc =
38 typename std::allocator_traits<Alloc>::template rebind_alloc<ScalarProd>;
39 auto sp = std::allocate_shared<ScalarProd>(ScalarProdAlloc(alloc));
40 return registry.create(type_name, op, sp, prec, config, alloc);
41 } else if constexpr (Concept::AssembledLinearOperator<Op>) {
43 return registry.create(type_name, op, config, alloc);
44 }
45 }
46
47 // not found: diagnose error
49 std::set<std::string> direct_keys;
50 if constexpr (Concept::AssembledLinearOperator<Op>)
52
54 "The key '{}' is not a known inverse operator type.\nPossible iterative "
55 "solvers: {}\nPossible direct solvers: {}",
59}
60
61} // namespace Impl
62
63template<Concept::LinearOperator Op, class Alloc = std::allocator<Op>>
64std::shared_ptr<InverseOperator<typename Op::domain_type, typename Op::range_type>>
65makeInverseOperator(const std::shared_ptr<Op>& op,
66 const ParameterTree& config,
67 const Alloc& alloc = Alloc())
68{
69 return Impl::makeInverseOperator(op, config, alloc);
70}
71
72} // Dune::Copasi::ISTL
73
74#endif // DUNE_COPASI_SOLVER_ISTL_FACTORY_INVERSE_HH
#define DUNE_COPASI_DEFAULT_LINEAR_INVERSE_OPERATOR
Definition inverse.hh:14
Definition block_jacobi.hh:14
std::shared_ptr< InverseOperator< typename Op::domain_type, typename Op::range_type > > makeInverseOperator(const std::shared_ptr< Op > &op, const ParameterTree &config, const Alloc &alloc=Alloc())
Definition inverse.hh:65
std::shared_ptr< Preconditioner< typename Op::domain_type, typename Op::range_type > > makePreconditioner(const std::shared_ptr< Op > &op, const ParameterTree &config, const Alloc &alloc=Alloc())
Definition preconditioner.hh:121
auto format_exception(Exception &&e, fmt::format_string< Args... > format, Args &&... args)
Definition exceptions.hh:23
constexpr bool is_bitflags_v
Alias for Bitflag indicator.
Definition bit_flags.hh:24