Loading [MathJax]/extensions/tex2jax.js

DUNE-FUNCTIONS (unstable)

istlvectorfactory.hh
1// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3// SPDX-FileCopyrightText: Copyright © DUNE Project contributors, see file AUTHORS.md
4// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception OR LGPL-3.0-or-later
5#ifndef DUNE_FUNCTIONS_BACKENDS_ISTL_VECTORFACTORY_HH
6#define DUNE_FUNCTIONS_BACKENDS_ISTL_VECTORFACTORY_HH
7
8#include <type_traits>
9
10#include <dune/common/exceptions.hh>
11#include <dune/common/fvector.hh>
12#include <dune/common/indices.hh>
13
14#include <dune/functions/functionspacebases/containerdescriptors.hh>
15
16#include <dune/istl/bvector.hh>
17#include <dune/istl/multitypeblockvector.hh>
18
19namespace Dune::Functions {
20namespace ContainerDescriptors {
21namespace Impl {
22
23template<class T>
24struct ISTLVectorFactory
25{
26 void operator() (const Unknown& tree) const
27 {
28 DUNE_THROW(Dune::NotImplemented, "Cannot create a vector. The container descriptor is unknown.");
29 }
30
31 template<class... V>
32 auto operator() (const Tuple<V...>& tree) const
33 {
34 return unpackIntegerSequence([&](auto... ii) {
35 return Dune::MultiTypeBlockVector<decltype((*this)(tree[ii]))...>{(*this)(tree[ii])...};
36 }, std::make_index_sequence<sizeof...(V)>());
37 }
38
39 template<class V, std::size_t n>
40 auto operator() (const Array<V,n>& tree) const
41 {
42 return unpackIntegerSequence([&](auto... ii) {
43 return Dune::BlockVector{(*this)(tree[ii])...};
44 }, std::make_index_sequence<n>());
45 }
46
47 template<class V>
48 auto operator() (const Vector<V>& tree) const
49 {
50 using W = decltype((*this)(tree[0]));
51 Dune::BlockVector<W> result(tree.size());
52 for (std::size_t i = 0; i < tree.size(); ++i)
53 result[i] = (*this)(tree[i]);
54 return result;
55 }
56
57 template<class V, std::size_t n>
58 auto operator() (const UniformArray<V,n>& tree) const
59 {
60 auto node = (*this)(tree[0]);
61 return unpackIntegerSequence([&](auto... ii) {
62 return Dune::BlockVector{((void)(ii),node)...};
63 }, std::make_index_sequence<n>());
64 }
65
66 template<class V>
67 auto operator() (const UniformVector<V>& tree) const
68 {
69 auto node = (*this)(tree[0]);
70 using W = decltype(node);
71 Dune::BlockVector<W> result(tree.size());
72 for (std::size_t i = 0; i < tree.size(); ++i)
73 result[i] = node;
74 return result;
75 }
76
77 // scalar types
78
79 auto operator() (const Value& tree) const
80 {
81 return T(0);
82 }
83
84 // flat vectors
85
86 template<std::size_t n>
87 auto operator() (const UniformArray<Value,n>& tree) const
88 {
89 return Dune::FieldVector<T,n>(0);
90 }
91
92 auto operator() (const UniformVector<Value>& tree) const
93 {
94 return Dune::BlockVector<T>(tree.size());
95 }
96
97 // block vectors
98
99 template<std::size_t n>
100 auto operator() (const UniformVector<UniformArray<Value,n>>& tree) const
101 {
102 return Dune::BlockVector<Dune::FieldVector<T,n>>(tree.size());
103 }
104
105 template<std::size_t n>
106 auto operator() (const Vector<UniformArray<Value,n>>& tree) const
107 {
108 return Dune::BlockVector<Dune::FieldVector<T,n>>(tree.size());
109 }
110
111 template<std::size_t n, std::size_t m>
112 auto operator() (const Array<UniformArray<Value,n>,m>& tree) const
113 {
114 return Dune::BlockVector<Dune::FieldVector<T,n>>(m);
115 }
116};
117
118} // end namespace Impl
119} // end namespace ContainerDescriptors
120
121
130template<class T = double, class ContainerDescriptor>
131auto makeISTLVector (const ContainerDescriptor& tree)
132{
133 auto factory = ContainerDescriptors::Impl::ISTLVectorFactory<T>{};
134 return factory(tree);
135}
136
137} // end namespace Dune::Functions
138
139#endif // DUNE_FUNCTIONS_BACKENDS_ISTL_VECTORFACTORY_HH
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden & Uni Heidelberg  |  generated with Hugo v0.111.3 (Mar 16, 23:47, 2025)