3 #ifndef DUNE_ISTL_MULTITYPEBLOCKVECTOR_HH
4 #define DUNE_ISTL_MULTITYPEBLOCKVECTOR_HH
7 #ifdef HAVE_BOOST_FUSION
12 #include <dune/common/dotproduct.hh>
13 #include <dune/common/ftraits.hh>
17 #include <boost/fusion/sequence.hpp>
18 #include <boost/fusion/container.hpp>
19 #include <boost/fusion/iterator.hpp>
20 #include <boost/typeof/typeof.hpp>
21 #include <boost/fusion/algorithm.hpp>
23 namespace mpl=boost::mpl;
24 namespace fusion=boost::fusion;
28 template<
typename T1,
typename T2=fusion::void_,
typename T3=fusion::void_,
typename T4=fusion::void_,
29 typename T5=fusion::void_,
typename T6=fusion::void_,
typename T7=fusion::void_,
30 typename T8=fusion::void_,
typename T9=fusion::void_>
31 class MultiTypeBlockVector;
61 template<
int current_element,
int remaining_elements,
typename TVec>
62 class MultiTypeBlockVector_Print {
67 static void print(
const TVec& v) {
68 std::cout <<
"\t(" << current_element <<
"):\n" << fusion::at_c<current_element>(v) <<
"\n";
69 MultiTypeBlockVector_Print<current_element+1,remaining_elements-1,TVec>::print(v);
72 template<
int current_element,
typename TVec>
73 class MultiTypeBlockVector_Print<current_element,0,TVec> {
75 static void print(
const TVec& v) {std::cout <<
"\n";}
91 template<
int count,
typename T1,
typename T2>
92 class MultiTypeBlockVector_Ident {
99 static void equalize(T1& a,
const T2& b) {
100 fusion::at_c<
count-1>(a) = b;
101 MultiTypeBlockVector_Ident<count-1,T1,T2>::equalize(a,b);
104 template<
typename T1,
typename T2>
105 class MultiTypeBlockVector_Ident<0,T1,T2> {
public:
static void equalize (T1& a,
const T2& b) {} };
117 template<
int count,
typename T>
118 class MultiTypeBlockVector_Add {
124 static void add (T& a,
const T& b) {
125 fusion::at_c<(
count-1)>(a) += fusion::at_c<(
count-1)>(b);
126 MultiTypeBlockVector_Add<count-1,T>::add(a,b);
132 static void sub (T& a,
const T& b) {
133 fusion::at_c<(
count-1)>(a) -= fusion::at_c<(
count-1)>(b);
134 MultiTypeBlockVector_Add<count-1,T>::sub(a,b);
138 class MultiTypeBlockVector_Add<0,T> {
public:
static void add (T& a,
const T& b) {}
static void sub (T& a,
const T& b) {} };
147 template<
int count,
typename TVec,
typename Ta>
148 class MultiTypeBlockVector_AXPY {
154 static void axpy(TVec& x,
const Ta& a,
const TVec& y) {
155 fusion::at_c<(
count-1)>(x).axpy(a,fusion::at_c<(
count-1)>(y));
156 MultiTypeBlockVector_AXPY<count-1,TVec,Ta>::axpy(x,a,y);
159 template<
typename TVec,
typename Ta>
160 class MultiTypeBlockVector_AXPY<0,TVec,Ta> {
public:
static void axpy (TVec& x,
const Ta& a,
const TVec& y) {} };
168 template<
int count,
typename TVec,
typename Ta>
169 class MultiTypeBlockVector_Mulscal {
175 static void mul(TVec& x,
const Ta& a) {
176 fusion::at_c<(
count-1)>(x) *= a;
177 MultiTypeBlockVector_Mulscal<count-1,TVec,Ta>::mul(x,a);
180 template<
typename TVec,
typename Ta>
181 class MultiTypeBlockVector_Mulscal<0,TVec,Ta> {
public:
static void mul (TVec& x,
const Ta& a) {} };
193 template<
int count,
typename TVec>
194 class MultiTypeBlockVector_Mul {
196 static typename TVec::field_type mul(
const TVec& x,
const TVec& y) {
return (fusion::at_c<count-1>(x) * fusion::at_c<count-1>(y)) + MultiTypeBlockVector_Mul<count-1,TVec>::mul(x,y); }
197 static typename TVec::field_type dot(
const TVec& x,
const TVec& y) {
return (Dune::dot(fusion::at_c<count-1>(x),fusion::at_c<count-1>(y))) + MultiTypeBlockVector_Mul<count-1,TVec>::dot(x,y); }
199 template<
typename TVec>
200 class MultiTypeBlockVector_Mul<0,TVec> {
202 static typename TVec::field_type mul(
const TVec& x,
const TVec& y) {
return 0;}
203 static typename TVec::field_type dot(
const TVec& x,
const TVec& y) {
return 0;}
216 template<
int count,
typename T>
217 class MultiTypeBlockVector_Norm {
219 typedef typename T::field_type field_type;
220 typedef typename FieldTraits<field_type>::real_type real_type;
225 static real_type result (
const T& a) {
226 return fusion::at_c<
count-1>(a).two_norm2() + MultiTypeBlockVector_Norm<count-1,T>::result(a);
231 class MultiTypeBlockVector_Norm<0,T> {
233 typedef typename T::field_type field_type;
234 typedef typename FieldTraits<field_type>::real_type real_type;
235 static real_type result (
const T& a) {
return 0.0;}
246 template<
typename T1,
typename T2,
typename T3,
typename T4,
247 typename T5,
typename T6,
typename T7,
typename T8,
typename T9>
248 class MultiTypeBlockVector :
public fusion::vector<T1, T2, T3, T4, T5, T6, T7, T8, T9> {
255 typedef MultiTypeBlockVector<T1, T2, T3, T4, T5, T6, T7, T8, T9> type;
257 typedef typename T1::field_type field_type;
262 int count() {
return mpl::size<type>::value;}
282 template< std::
size_t index >
284 operator[] (
const std::integral_constant< std::size_t, index > indexVariable )
285 -> decltype(fusion::at_c<index>(*
this))
287 DUNE_UNUSED_PARAMETER(indexVariable);
288 return fusion::at_c<index>(*this);
296 template< std::
size_t index >
298 operator[] (
const std::integral_constant< std::size_t, index > indexVariable )
const
299 -> decltype(fusion::at_c<index>(*
this))
301 DUNE_UNUSED_PARAMETER(indexVariable);
302 return fusion::at_c<index>(*this);
309 void operator= (
const T& newval) {MultiTypeBlockVector_Ident<mpl::size<type>::value,type,T>::equalize(*
this, newval); }
314 void operator+= (
const type& newv) {MultiTypeBlockVector_Add<mpl::size<type>::value,type>::add(*
this,newv);}
319 void operator-= (
const type& newv) {MultiTypeBlockVector_Add<mpl::size<type>::value,type>::sub(*
this,newv);}
321 void operator*= (
const int& w) {MultiTypeBlockVector_Mulscal<mpl::size<type>::value,type,
const int>::mul(*
this,w);}
322 void operator*= (
const float& w) {MultiTypeBlockVector_Mulscal<mpl::size<type>::value,type,
const float>::mul(*
this,w);}
323 void operator*= (
const double& w) {MultiTypeBlockVector_Mulscal<mpl::size<type>::value,type,
const double>::mul(*
this,w);}
325 field_type operator* (
const type& newv)
const {
return MultiTypeBlockVector_Mul<mpl::size<type>::value,type>::mul(*
this,newv);}
326 field_type dot (
const type& newv)
const {
return MultiTypeBlockVector_Mul<mpl::size<type>::value,type>::dot(*
this,newv);}
331 typename FieldTraits<field_type>::real_type two_norm2()
const {
return MultiTypeBlockVector_Norm<mpl::size<type>::value,type>::result(*
this);}
336 typename FieldTraits<field_type>::real_type two_norm()
const {
return sqrt(this->two_norm2());}
341 template<
typename Ta>
342 void axpy (
const Ta& a,
const type& y) {
343 MultiTypeBlockVector_AXPY<mpl::size<type>::value,type,Ta>::axpy(*
this,a,y);
355 template<
typename T1,
typename T2,
typename T3,
typename T4,
typename T5,
typename T6,
typename T7,
typename T8,
typename T9>
356 std::ostream& operator<< (std::ostream& s, const MultiTypeBlockVector<T1,T2,T3,T4,T5,T6,T7,T8,T9>& v) {
357 MultiTypeBlockVector_Print<0,mpl::size<MultiTypeBlockVector<T1,T2,T3,T4,T5,T6,T7,T8,T9> >::value,MultiTypeBlockVector<T1,T2,T3,T4,T5,T6,T7,T8,T9> >::print(v);
365 #endif // end HAVE_BOOST_FUSION
366 #endif // end HAVE_DUNE_BOOST
std::size_t count
Definition: matrixmatrix.hh:215
Simple iterative methods like Jacobi, Gauss-Seidel, SOR, SSOR, etc. in a generic way.