Skip to content
Snippets Groups Projects
Commit 1462583d authored by Oliver Sander's avatar Oliver Sander
Browse files

Reimplement MultiTypeBlockVector using std::tuple instead of boost::fusion

This removes the dependency on boost::fusion.  Using a bit of variadic template
magic makes the code a tiny bit simpler, too.

Note, though, that this change is not backward-compatible.  If you have previously
used

boost::fusion::at_c<i>(v)

to get the i-th entry of the vector v, you now have to change this to

std::get<i>(v)

Alternatively, you can start using operator[], and write

v[std::integral_constant<int,i>()]
parent 3acd7af5
No related branches found
No related tags found
No related merge requests found
......@@ -378,12 +378,12 @@ namespace Dune {
template<typename T11, typename T12, typename T13, typename T14,
typename T15, typename T16, typename T17, typename T18,
typename T19, typename T21, typename T22, typename T23,
typename T24, typename T25, typename T26, typename T27,
typename T28, typename T29, class K>
typename T19,
typename... MultiTypeVectorArgs,
class K>
static void dbgs (const MultiTypeBlockMatrix<T11,T12,T13,T14,T15,T16,T17,T18,T19>& A,
MultiTypeBlockVector<T21,T22,T23,T24,T25,T26,T27,T28,T29>& x,
const MultiTypeBlockVector<T21,T22,T23,T24,T25,T26,T27,T28,T29>& b,
MultiTypeBlockVector<MultiTypeVectorArgs...>& x,
const MultiTypeBlockVector<MultiTypeVectorArgs...>& b,
const K& w)
{
const int rowcount =
......@@ -426,12 +426,12 @@ namespace Dune {
template<typename T11, typename T12, typename T13, typename T14,
typename T15, typename T16, typename T17, typename T18,
typename T19, typename T21, typename T22, typename T23,
typename T24, typename T25, typename T26, typename T27,
typename T28, typename T29, class K>
typename T19,
typename... MultiTypeVectorArgs,
class K>
static void bsorf (const MultiTypeBlockMatrix<T11,T12,T13,T14,T15,T16,T17,T18,T19>& A,
MultiTypeBlockVector<T21,T22,T23,T24,T25,T26,T27,T28,T29>& x,
const MultiTypeBlockVector<T21,T22,T23,T24,T25,T26,T27,T28,T29>& b,
MultiTypeBlockVector<MultiTypeVectorArgs...>& x,
const MultiTypeBlockVector<MultiTypeVectorArgs...>& b,
const K& w)
{
const int rowcount =
......@@ -476,12 +476,12 @@ namespace Dune {
template<typename T11, typename T12, typename T13, typename T14,
typename T15, typename T16, typename T17, typename T18,
typename T19, typename T21, typename T22, typename T23,
typename T24, typename T25, typename T26, typename T27,
typename T28, typename T29, class K>
typename T19,
typename... MultiTypeVectorArgs,
class K>
static void bsorb (const MultiTypeBlockMatrix<T11,T12,T13,T14,T15,T16,T17,T18,T19>& A,
MultiTypeBlockVector<T21,T22,T23,T24,T25,T26,T27,T28,T29>& x,
const MultiTypeBlockVector<T21,T22,T23,T24,T25,T26,T27,T28,T29>& b,
MultiTypeBlockVector<MultiTypeVectorArgs...>& x,
const MultiTypeBlockVector<MultiTypeVectorArgs...>& b,
const K& w)
{
const int rowcount =
......@@ -523,15 +523,14 @@ namespace Dune {
#if HAVE_DUNE_BOOST
#ifdef HAVE_BOOST_FUSION
template<typename T11, typename T12, typename T13, typename T14,
typename T15, typename T16, typename T17, typename T18,
typename T19, typename T21, typename T22, typename T23,
typename T24, typename T25, typename T26, typename T27,
typename T28, typename T29, class K>
typename T19,
typename... MultiTypeVectorArgs,
class K>
static void dbjac (const MultiTypeBlockMatrix<T11,T12,T13,T14,T15,T16,T17,T18,T19>& A,
MultiTypeBlockVector<T21,T22,T23,T24,T25,T26,T27,T28,T29>& x,
const MultiTypeBlockVector<T21,T22,T23,T24,T25,T26,T27,T28,T29>& b,
MultiTypeBlockVector<MultiTypeVectorArgs...>& x,
const MultiTypeBlockVector<MultiTypeVectorArgs...>& b,
const K& w)
{
const int rowcount =
......
......@@ -135,14 +135,14 @@ namespace Dune {
/** \brief y += A x
*/
static void umv(TVecY& y, const TMatrix& A, const TVecX& x) {
fusion::at_c<ccol>( fusion::at_c<crow>(A) ).umv( fusion::at_c<ccol>(x), fusion::at_c<crow>(y) );
std::get<ccol>( fusion::at_c<crow>(A) ).umv( std::get<ccol>(x), std::get<crow>(y) );
MultiTypeBlockMatrix_VectMul<crow,remain_rows,ccol+1,remain_cols-1,TVecY,TMatrix,TVecX>::umv(y, A, x);
}
/** \brief y -= A x
*/
static void mmv(TVecY& y, const TMatrix& A, const TVecX& x) {
fusion::at_c<ccol>( fusion::at_c<crow>(A) ).mmv( fusion::at_c<ccol>(x), fusion::at_c<crow>(y) );
std::get<ccol>( fusion::at_c<crow>(A) ).mmv( std::get<ccol>(x), std::get<crow>(y) );
MultiTypeBlockMatrix_VectMul<crow,remain_rows,ccol+1,remain_cols-1,TVecY,TMatrix,TVecX>::mmv(y, A, x);
}
......@@ -151,7 +151,7 @@ namespace Dune {
*/
template<typename AlphaType>
static void usmv(const AlphaType& alpha, TVecY& y, const TMatrix& A, const TVecX& x) {
fusion::at_c<ccol>( fusion::at_c<crow>(A) ).usmv(alpha, fusion::at_c<ccol>(x), fusion::at_c<crow>(y) );
std::get<ccol>( fusion::at_c<crow>(A) ).usmv(alpha, std::get<ccol>(x), std::get<crow>(y) );
MultiTypeBlockMatrix_VectMul<crow,remain_rows,ccol+1,remain_cols-1,TVecY,TMatrix,TVecX>::usmv(alpha,y, A, x);
}
......@@ -224,7 +224,7 @@ namespace Dune {
*/
typedef MultiTypeBlockMatrix<T1, T2, T3, T4, T5, T6, T7, T8, T9> type;
typedef typename mpl::at_c<T1,0>::type field_type;
typedef typename T1::field_type field_type;
/**
* assignment operator
......@@ -236,7 +236,7 @@ namespace Dune {
*/
template<typename X, typename Y>
void mv (const X& x, Y& y) const {
static_assert(x.size() == mpl::size<T1>::value, "length of x does not match row length");
static_assert(x.size() == T1::size(), "length of x does not match row length");
static_assert(y.size() == mpl::size<type>::value, "length of y does not match row count");
y = 0; //reset y (for mv uses umv)
......@@ -247,7 +247,7 @@ namespace Dune {
*/
template<typename X, typename Y>
void umv (const X& x, Y& y) const {
static_assert(x.size() == mpl::size<T1>::value, "length of x does not match row length");
static_assert(x.size() == T1::size(), "length of x does not match row length");
static_assert(y.size() == mpl::size<type>::value, "length of y does not match row count");
MultiTypeBlockMatrix_VectMul<0,mpl::size<type>::value,0,T1::size(),Y,type,X>::umv(y, *this, x); //iterate over all matrix elements
......@@ -257,7 +257,7 @@ namespace Dune {
*/
template<typename X, typename Y>
void mmv (const X& x, Y& y) const {
static_assert(x.size() == mpl::size<T1>::value, "length of x does not match row length");
static_assert(x.size() == T1::size(), "length of x does not match row length");
static_assert(y.size() == mpl::size<type>::value, "length of y does not match row count");
MultiTypeBlockMatrix_VectMul<0,mpl::size<type>::value,0,T1::size(),Y,type,X>::mmv(y, *this, x); //iterate over all matrix elements
......@@ -267,7 +267,7 @@ namespace Dune {
*/
template<typename AlphaType, typename X, typename Y>
void usmv (const AlphaType& alpha, const X& x, Y& y) const {
static_assert(x.size() == mpl::size<T1>::value, "length of x does not match row length");
static_assert(x.size() == T1::size(), "length of x does not match row length");
static_assert(y.size() == mpl::size<type>::value, "length of y does not match row count");
MultiTypeBlockMatrix_VectMul<0,mpl::size<type>::value,0,T1::size(),Y,type,X>::usmv(alpha,y, *this, x); //iterate over all matrix elements
......@@ -321,7 +321,7 @@ namespace Dune {
*/
template <typename Trhs, typename TVector, typename TMatrix, typename K>
static void calc_rhs(const TMatrix& A, TVector& x, TVector& v, Trhs& b, const K& w) {
fusion::at_c<ccol>( fusion::at_c<crow>(A) ).mmv( fusion::at_c<ccol>(x), b );
std::get<ccol>( fusion::at_c<crow>(A) ).mmv( std::get<ccol>(x), b );
MultiTypeBlockMatrix_Solver_Col<I, crow, ccol+1, remain_col-1>::calc_rhs(A,x,v,b,w); //next column element
}
......@@ -358,12 +358,11 @@ namespace Dune {
}
template <typename TVector, typename TMatrix, typename K>
static void dbgs(const TMatrix& A, TVector& x, TVector& v, const TVector& b, const K& w) {
typename mpl::at_c<TVector,crow>::type rhs;
rhs = fusion::at_c<crow> (b);
auto rhs = std::get<crow> (b);
MultiTypeBlockMatrix_Solver_Col<I,crow,0, mpl::at_c<TMatrix,crow>::type::size()>::calc_rhs(A,x,v,rhs,w); // calculate right side of equation
//solve on blocklevel I-1
algmeta_itsteps<I-1>::dbgs(fusion::at_c<crow>( fusion::at_c<crow>(A)), fusion::at_c<crow>(x),rhs,w);
algmeta_itsteps<I-1>::dbgs(std::get<crow>( fusion::at_c<crow>(A)), std::get<crow>(x),rhs,w);
MultiTypeBlockMatrix_Solver<I,crow+1,remain_row-1>::dbgs(A,x,v,b,w); //next row
}
......@@ -381,13 +380,12 @@ namespace Dune {
}
template <typename TVector, typename TMatrix, typename K> //recursion over all matrix rows (A)
static void bsorf(const TMatrix& A, TVector& x, TVector& v, const TVector& b, const K& w) {
typename mpl::at_c<TVector,crow>::type rhs;
rhs = fusion::at_c<crow> (b);
auto rhs = std::get<crow> (b);
MultiTypeBlockMatrix_Solver_Col<I,crow,0, mpl::at_c<TMatrix,crow>::type::size()>::calc_rhs(A,x,v,rhs,w); // calculate right side of equation
//solve on blocklevel I-1
algmeta_itsteps<I-1>::bsorf(fusion::at_c<crow>( fusion::at_c<crow>(A)), fusion::at_c<crow>(v),rhs,w);
fusion::at_c<crow>(x).axpy(w,fusion::at_c<crow>(v));
algmeta_itsteps<I-1>::bsorf(std::get<crow>( fusion::at_c<crow>(A)), std::get<crow>(v),rhs,w);
std::get<crow>(x).axpy(w,std::get<crow>(v));
MultiTypeBlockMatrix_Solver<I,crow+1,remain_row-1>::bsorf(A,x,v,b,w); //next row
}
......@@ -403,13 +401,12 @@ namespace Dune {
}
template <typename TVector, typename TMatrix, typename K> //recursion over all matrix rows (A)
static void bsorb(const TMatrix& A, TVector& x, TVector& v, const TVector& b, const K& w) {
typename mpl::at_c<TVector,crow>::type rhs;
rhs = fusion::at_c<crow> (b);
auto rhs = std::get<crow> (b);
MultiTypeBlockMatrix_Solver_Col<I,crow,0, mpl::at_c<TMatrix,crow>::type::size()>::calc_rhs(A,x,v,rhs,w); // calculate right side of equation
//solve on blocklevel I-1
algmeta_itsteps<I-1>::bsorb(fusion::at_c<crow>( fusion::at_c<crow>(A)), fusion::at_c<crow>(v),rhs,w);
fusion::at_c<crow>(x).axpy(w,fusion::at_c<crow>(v));
algmeta_itsteps<I-1>::bsorb(std::get<crow>( fusion::at_c<crow>(A)), std::get<crow>(v),rhs,w);
std::get<crow>(x).axpy(w,std::get<crow>(v));
MultiTypeBlockMatrix_Solver<I,crow-1,remain_row-1>::bsorb(A,x,v,b,w); //next row
}
......@@ -426,12 +423,11 @@ namespace Dune {
}
template <typename TVector, typename TMatrix, typename K>
static void dbjac(const TMatrix& A, TVector& x, TVector& v, const TVector& b, const K& w) {
typename mpl::at_c<TVector,crow>::type rhs;
rhs = fusion::at_c<crow> (b);
auto rhs = std::get<crow> (b);
MultiTypeBlockMatrix_Solver_Col<I,crow,0, mpl::at_c<TMatrix,crow>::type::size()>::calc_rhs(A,x,v,rhs,w); // calculate right side of equation
//solve on blocklevel I-1
algmeta_itsteps<I-1>::dbjac(fusion::at_c<crow>( fusion::at_c<crow>(A)), fusion::at_c<crow>(v),rhs,w);
algmeta_itsteps<I-1>::dbjac(std::get<crow>( fusion::at_c<crow>(A)), std::get<crow>(v),rhs,w);
MultiTypeBlockMatrix_Solver<I,crow+1,remain_row-1>::dbjac(A,x,v,b,w); //next row
}
......
......@@ -3,31 +3,18 @@
#ifndef DUNE_ISTL_MULTITYPEBLOCKVECTOR_HH
#define DUNE_ISTL_MULTITYPEBLOCKVECTOR_HH
#if HAVE_DUNE_BOOST
#ifdef HAVE_BOOST_FUSION
#include <cmath>
#include <iostream>
#include <tuple>
#include <dune/common/dotproduct.hh>
#include <dune/common/ftraits.hh>
#include "istlexception.hh"
#include <boost/fusion/sequence.hpp>
#include <boost/fusion/container.hpp>
#include <boost/fusion/iterator.hpp>
#include <boost/typeof/typeof.hpp>
#include <boost/fusion/algorithm.hpp>
namespace mpl=boost::mpl;
namespace fusion=boost::fusion;
// forward declaration
namespace Dune {
template<typename T1, typename T2=fusion::void_, typename T3=fusion::void_, typename T4=fusion::void_,
typename T5=fusion::void_, typename T6=fusion::void_, typename T7=fusion::void_,
typename T8=fusion::void_, typename T9=fusion::void_>
template < typename... Args >
class MultiTypeBlockVector;
}
......@@ -65,7 +52,7 @@ namespace Dune {
* print out the current vector element and all following
*/
static void print(const TVec& v) {
std::cout << "\t(" << current_element << "):\n" << fusion::at_c<current_element>(v) << "\n";
std::cout << "\t(" << current_element << "):\n" << std::get<current_element>(v) << "\n";
MultiTypeBlockVector_Print<current_element+1,remaining_elements-1,TVec>::print(v); //next element
}
};
......@@ -97,7 +84,7 @@ namespace Dune {
* note: each MultiTypeBlockVector element has to provide the = operator with type T2
*/
static void equalize(T1& a, const T2& b) {
fusion::at_c<count-1>(a) = b; //equalize current elements
std::get<count-1>(a) = b; //equalize current elements
MultiTypeBlockVector_Ident<count-1,T1,T2>::equalize(a,b); //next elements
}
};
......@@ -122,7 +109,7 @@ namespace Dune {
* add vector to vector
*/
static void add (T& a, const T& b) { //add vector elements
fusion::at_c<(count-1)>(a) += fusion::at_c<(count-1)>(b);
std::get<(count-1)>(a) += std::get<(count-1)>(b);
MultiTypeBlockVector_Add<count-1,T>::add(a,b);
}
......@@ -130,7 +117,7 @@ namespace Dune {
* Subtract vector from vector
*/
static void sub (T& a, const T& b) { //sub vector elements
fusion::at_c<(count-1)>(a) -= fusion::at_c<(count-1)>(b);
std::get<(count-1)>(a) -= std::get<(count-1)>(b);
MultiTypeBlockVector_Add<count-1,T>::sub(a,b);
}
};
......@@ -152,7 +139,7 @@ namespace Dune {
* calculates x += a * y
*/
static void axpy(TVec& x, const Ta& a, const TVec& y) {
fusion::at_c<(count-1)>(x).axpy(a,fusion::at_c<(count-1)>(y));
std::get<(count-1)>(x).axpy(a,std::get<(count-1)>(y));
MultiTypeBlockVector_AXPY<count-1,TVec,Ta>::axpy(x,a,y);
}
};
......@@ -172,7 +159,7 @@ namespace Dune {
* calculates x *= a
*/
static void mul(TVec& x, const Ta& a) {
fusion::at_c<(count-1)>(x) *= a;
std::get<(count-1)>(x) *= a;
MultiTypeBlockVector_Mulscal<count-1,TVec,Ta>::mul(x,a);
}
};
......@@ -191,9 +178,17 @@ namespace Dune {
template<int count, typename TVec>
class MultiTypeBlockVector_Mul {
public:
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); }
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); }
static typename TVec::field_type mul(const TVec& x, const TVec& y)
{
return (std::get<count-1>(x) * std::get<count-1>(y)) + MultiTypeBlockVector_Mul<count-1,TVec>::mul(x,y);
}
static typename TVec::field_type dot(const TVec& x, const TVec& y)
{
return Dune::dot(std::get<count-1>(x),std::get<count-1>(y)) + MultiTypeBlockVector_Mul<count-1,TVec>::dot(x,y);
}
};
template<typename TVec>
class MultiTypeBlockVector_Mul<0,TVec> {
public:
......@@ -220,7 +215,7 @@ namespace Dune {
* sum up all elements' 2-norms
*/
static real_type result (const T& a) { //result = sum of all elements' 2-norms
return fusion::at_c<count-1>(a).two_norm2() + MultiTypeBlockVector_Norm<count-1,T>::result(a);
return std::get<count-1>(a).two_norm2() + MultiTypeBlockVector_Norm<count-1,T>::result(a);
}
};
......@@ -236,33 +231,42 @@ namespace Dune {
@brief A Vector class to support different block types
This vector class combines elements of different types known at compile-time.
You must add BOOST_CPPFLAGS and BOOT_LDFLAGS to the CPPFLAGS and LDFLAGS during
compilation, respectively, to use this class
*/
template<typename T1, typename T2, typename T3, typename T4,
typename T5, typename T6, typename T7, typename T8, typename T9>
class MultiTypeBlockVector : public fusion::vector<T1, T2, T3, T4, T5, T6, T7, T8, T9> {
template < typename... Args >
class MultiTypeBlockVector
: public std::tuple<Args...>
{
/** \brief Helper type */
typedef std::tuple<Args...> tupleType;
public:
/**
* own class' type
*/
typedef MultiTypeBlockVector<T1, T2, T3, T4, T5, T6, T7, T8, T9> type;
typedef MultiTypeBlockVector<Args...> type;
typedef typename T1::field_type field_type;
/** \brief The type used for scalars
*
* The current code hardwires it to 'double', which is far from nice.
* On the other hand, it is not clear what the correct type is. If the MultiTypeBlockVector class
* is instantiated with several vectors of different field_types, what should the resulting
* field_type be?
*/
typedef double field_type;
/** \brief Return the number of vector entries */
static DUNE_CONSTEXPR std::size_t size()
{
return mpl::size<type>::value;
return sizeof...(Args);
}
/**
* number of elements
*/
int count() {return mpl::size<type>::value;}
int count()
{
return sizeof...(Args);
}
/** \brief Random-access operator
*
......@@ -283,11 +287,11 @@ namespace Dune {
* \endcode
*/
template< int index >
typename mpl::at_c<type,index>::type&
typename std::tuple_element<index,tupleType>::type&
operator[] ( const std::integral_constant< int, index > indexVariable )
{
DUNE_UNUSED_PARAMETER(indexVariable);
return fusion::at_c<index>(*this);
return std::get<index>(*this);
}
/** \brief Const random-access operator
......@@ -296,38 +300,38 @@ namespace Dune {
* explanation of how to use it.
*/
template< int index >
const typename mpl::at_c<type,index>::type&
const typename std::tuple_element<index,tupleType>::type&
operator[] ( const std::integral_constant< int, index > indexVariable ) const
{
DUNE_UNUSED_PARAMETER(indexVariable);
return fusion::at_c<index>(*this);
return std::get<index>(*this);
}
/** \brief Assignment operator
*/
template<typename T>
void operator= (const T& newval) {MultiTypeBlockVector_Ident<mpl::size<type>::value,type,T>::equalize(*this, newval); }
void operator= (const T& newval) {MultiTypeBlockVector_Ident<sizeof...(Args),type,T>::equalize(*this, newval); }
/**
* operator for MultiTypeBlockVector += MultiTypeBlockVector operations
*/
void operator+= (const type& newv) {MultiTypeBlockVector_Add<mpl::size<type>::value,type>::add(*this,newv);}
void operator+= (const type& newv) {MultiTypeBlockVector_Add<sizeof...(Args),type>::add(*this,newv);}
/**
* operator for MultiTypeBlockVector -= MultiTypeBlockVector operations
*/
void operator-= (const type& newv) {MultiTypeBlockVector_Add<mpl::size<type>::value,type>::sub(*this,newv);}
void operator-= (const type& newv) {MultiTypeBlockVector_Add<sizeof...(Args),type>::sub(*this,newv);}
void operator*= (const int& w) {MultiTypeBlockVector_Mulscal<mpl::size<type>::value,type,const int>::mul(*this,w);}
void operator*= (const float& w) {MultiTypeBlockVector_Mulscal<mpl::size<type>::value,type,const float>::mul(*this,w);}
void operator*= (const double& w) {MultiTypeBlockVector_Mulscal<mpl::size<type>::value,type,const double>::mul(*this,w);}
void operator*= (const int& w) {MultiTypeBlockVector_Mulscal<sizeof...(Args),type,const int>::mul(*this,w);}
void operator*= (const float& w) {MultiTypeBlockVector_Mulscal<sizeof...(Args),type,const float>::mul(*this,w);}
void operator*= (const double& w) {MultiTypeBlockVector_Mulscal<sizeof...(Args),type,const double>::mul(*this,w);}
field_type operator* (const type& newv) const {return MultiTypeBlockVector_Mul<mpl::size<type>::value,type>::mul(*this,newv);}
field_type dot (const type& newv) const {return MultiTypeBlockVector_Mul<mpl::size<type>::value,type>::dot(*this,newv);}
field_type operator* (const type& newv) const {return MultiTypeBlockVector_Mul<sizeof...(Args),type>::mul(*this,newv);}
field_type dot (const type& newv) const {return MultiTypeBlockVector_Mul<sizeof...(Args),type>::dot(*this,newv);}
/** \brief Compute the squared Euclidean norm
*/
typename FieldTraits<field_type>::real_type two_norm2() const {return MultiTypeBlockVector_Norm<mpl::size<type>::value,type>::result(*this);}
typename FieldTraits<field_type>::real_type two_norm2() const {return MultiTypeBlockVector_Norm<sizeof...(Args),type>::result(*this);}
/** \brief Compute the Euclidean norm
*/
......@@ -339,7 +343,7 @@ namespace Dune {
*/
template<typename Ta>
void axpy (const Ta& a, const type& y) {
MultiTypeBlockVector_AXPY<mpl::size<type>::value,type,Ta>::axpy(*this,a,y);
MultiTypeBlockVector_AXPY<sizeof...(Args),type,Ta>::axpy(*this,a,y);
}
};
......@@ -348,9 +352,9 @@ namespace Dune {
/** \brief Send MultiTypeBlockVector to an outstream
*/
template<typename T1, typename T2, typename T3, typename T4, typename T5, typename T6, typename T7, typename T8, typename T9>
std::ostream& operator<< (std::ostream& s, const MultiTypeBlockVector<T1,T2,T3,T4,T5,T6,T7,T8,T9>& v) {
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);
template <typename... Args>
std::ostream& operator<< (std::ostream& s, const MultiTypeBlockVector<Args...>& v) {
MultiTypeBlockVector_Print<0,sizeof...(Args),MultiTypeBlockVector<Args...> >::print(v);
return s;
}
......@@ -358,7 +362,4 @@ namespace Dune {
} // end namespace
#endif // end HAVE_BOOST_FUSION
#endif // end HAVE_DUNE_BOOST
#endif
......@@ -69,7 +69,6 @@ add_executable(mmtest mmtest.cc)
add_executable(multitypeblockmatrixtest "multitypeblockmatrixtest.cc")
add_dune_boost_flags(multitypeblockmatrixtest)
add_executable(multitypeblockvectortest "multitypeblockvectortest.cc")
add_dune_boost_flags(multitypeblockvectortest)
add_executable(mv "mv.cc")
add_executable(iotest "iotest.cc")
add_executable(inverseoperator2prectest "inverseoperator2prectest.cc")
......@@ -131,7 +130,7 @@ foreach(_exe ${ALLTESTS})
add_test(${_exe} ${_exe})
endforeach(_exe ${ALLTESTS})
set_tests_properties(multitypeblockmatrixtest multitypeblockvectortest PROPERTIES SKIP_RETURN_CODE 77)
set_tests_properties(multitypeblockmatrixtest PROPERTIES SKIP_RETURN_CODE 77)
exclude_from_headercheck(
complexdata.hh)
......@@ -41,37 +41,40 @@ int main(int argc, char** argv) try
typedef MultiTypeBlockVector<Matrix<FieldMatrix<double,3,3> >, Matrix<FieldMatrix<double,3,1> > > RowType0;
typedef MultiTypeBlockVector<Matrix<FieldMatrix<double,1,3> >, Matrix<FieldMatrix<double,1,1> > > RowType1;
std::integral_constant<int, 0> _0;
std::integral_constant<int, 1> _1;
MultiTypeBlockMatrix<RowType0,RowType1> multiMatrix;
fusion::at_c<0>( fusion::at_c<0>(multiMatrix)).setSize(3,3);
fusion::at_c<1>( fusion::at_c<0>(multiMatrix)).setSize(3,2);
fusion::at_c<0>( fusion::at_c<1>(multiMatrix)).setSize(2,3);
fusion::at_c<1>( fusion::at_c<1>(multiMatrix)).setSize(2,2);
fusion::at_c<0>(multiMatrix)[_0].setSize(3,3);
fusion::at_c<0>(multiMatrix)[_1].setSize(3,2);
fusion::at_c<1>(multiMatrix)[_0].setSize(2,3);
fusion::at_c<1>(multiMatrix)[_1].setSize(2,2);
// lazy solution: initialize the entire matrix with zeros
fusion::at_c<0>( fusion::at_c<0>(multiMatrix)) = 0;
fusion::at_c<1>( fusion::at_c<0>(multiMatrix)) = 0;
fusion::at_c<0>( fusion::at_c<1>(multiMatrix)) = 0;
fusion::at_c<1>( fusion::at_c<1>(multiMatrix)) = 0;
fusion::at_c<0>(multiMatrix)[_0] = 0;
fusion::at_c<0>(multiMatrix)[_1] = 0;
fusion::at_c<1>(multiMatrix)[_0] = 0;
fusion::at_c<1>(multiMatrix)[_1] = 0;
printmatrix(std::cout, fusion::at_c<0>( fusion::at_c<0>(multiMatrix)), "(0,0)", "--");
printmatrix(std::cout, fusion::at_c<1>( fusion::at_c<0>(multiMatrix)), "(0,1)", "--");
printmatrix(std::cout, fusion::at_c<0>( fusion::at_c<1>(multiMatrix)), "(1,0)", "--");
printmatrix(std::cout, fusion::at_c<1>( fusion::at_c<1>(multiMatrix)), "(1,1)", "--");
printmatrix(std::cout, fusion::at_c<0>(multiMatrix)[_0], "(0,0)", "--");
printmatrix(std::cout, fusion::at_c<0>(multiMatrix)[_1], "(0,1)", "--");
printmatrix(std::cout, fusion::at_c<1>(multiMatrix)[_0], "(1,0)", "--");
printmatrix(std::cout, fusion::at_c<1>(multiMatrix)[_1], "(1,1)", "--");
// set up a test vector
MultiTypeBlockVector<BlockVector<FieldVector<double,3> >, BlockVector<FieldVector<double,1> > > multiVector;
boost::fusion::at_c<0>(multiVector) = {{1,0,0},
{0,1,0},
{0,0,1}};
multiVector[_0] = {{1,0,0},
{0,1,0},
{0,0,1}};
boost::fusion::at_c<1>(multiVector) = {3.14, 42};
multiVector[_1] = {3.14, 42};
// Test matrix-vector products
MultiTypeBlockVector<BlockVector<FieldVector<double,3> >, BlockVector<FieldVector<double,1> > > result;
fusion::at_c<0>(result).resize(3);
fusion::at_c<1>(result).resize(2);
result[_0].resize(3);
result[_1].resize(2);
multiMatrix.mv(multiVector,result);
multiMatrix.umv(multiVector,result);
......@@ -100,10 +103,10 @@ int main(int argc, char** argv) try
typedef Dune::MultiTypeBlockVector<Dune::BlockVector<Dune::FieldVector<double,1> >,Dune::BlockVector<Dune::FieldVector<double,1> > > TestVector;
TestVector x, b;
fusion::at_c<0>(x).resize(Y1);
fusion::at_c<1>(x).resize(Y2);
fusion::at_c<0>(b).resize(X1);
fusion::at_c<1>(b).resize(X2);
x[_0].resize(Y1);
x[_1].resize(Y2);
b[_0].resize(X1);
b[_1].resize(X2);
x = 1; b = 1;
......@@ -171,10 +174,10 @@ int main(int argc, char** argv) try
typedef Dune::MultiTypeBlockVector<BCRSMat,BCRSMat> BCRS_Row;
typedef Dune::MultiTypeBlockMatrix<BCRS_Row,BCRS_Row> CM_BCRS;
CM_BCRS A;
fusion::at_c<0>(fusion::at_c<0>(A)) = A11;
fusion::at_c<1>(fusion::at_c<0>(A)) = A12;
fusion::at_c<0>(fusion::at_c<1>(A)) = A21;
fusion::at_c<1>(fusion::at_c<1>(A)) = A22;
fusion::at_c<0>(A)[_0] = A11;
fusion::at_c<0>(A)[_1] = A12;
fusion::at_c<1>(A)[_0] = A21;
fusion::at_c<1>(A)[_1] = A22;
x = 1;
b = 1;
......@@ -191,8 +194,8 @@ int main(int argc, char** argv) try
// Solve linear system
loop.apply(x,b,r);
printvector(std::cout,fusion::at_c<0>(x),"solution x1","entry",11,9,1);
printvector(std::cout,fusion::at_c<1>(x),"solution x2","entry",11,9,1);
printvector(std::cout,x[_0],"solution x1","entry",11,9,1);
printvector(std::cout,x[_1],"solution x2","entry",11,9,1);
return 0;
#endif
......
......@@ -21,9 +21,6 @@ using namespace Dune;
int main(int argc, char** argv) try
{
#if ! HAVE_DUNE_BOOST || ! defined HAVE_BOOST_FUSION
return 77;
#else
MultiTypeBlockVector<BlockVector<FieldVector<double,3> >, BlockVector<FieldVector<double,1> > > multiVector;
std::integral_constant<int, 0> _0;
......@@ -79,7 +76,6 @@ int main(int argc, char** argv) try
std::cout << multiVector.dot(multiVector2) << std::endl;
return 0;
#endif
}
catch (Dune::Exception& e)
{
......
......@@ -516,6 +516,10 @@ void test_MultiTypeBlockVector_MultiTypeBlockMatrix() {
typedef Dune::FieldMatrix<double,1,1> LittleBlock; //matrix block type
typedef Dune::BCRSMatrix<LittleBlock> BCRSMat; //matrix type
// Compile-time constants for the indices '0' and '1'
std::integral_constant<int, 0> _0;
std::integral_constant<int, 1> _1;
const int X1=3; //index bounds of all four matrices
const int X2=2;
const int Y1=3;
......@@ -528,10 +532,10 @@ void test_MultiTypeBlockVector_MultiTypeBlockMatrix() {
typedef Dune::MultiTypeBlockVector<Dune::BlockVector<Dune::FieldVector<double,1> >,Dune::BlockVector<Dune::FieldVector<double,1> > > TestVector;
TestVector x, b;
fusion::at_c<0>(x).resize(Y1);
fusion::at_c<1>(x).resize(Y2);
fusion::at_c<0>(b).resize(X1);
fusion::at_c<1>(b).resize(X2);
x[_0].resize(Y1);
x[_1].resize(Y2);
b[_0].resize(X1);
b[_1].resize(X2);
x = 1; b = 1;
......@@ -568,10 +572,10 @@ void test_MultiTypeBlockVector_MultiTypeBlockMatrix() {
typedef Dune::MultiTypeBlockVector<BCRSMat,BCRSMat> BCRS_Row;
typedef Dune::MultiTypeBlockMatrix<BCRS_Row,BCRS_Row> CM_BCRS;
CM_BCRS A;
fusion::at_c<0>(fusion::at_c<0>(A)) = A11;
fusion::at_c<1>(fusion::at_c<0>(A)) = A12;
fusion::at_c<0>(fusion::at_c<1>(A)) = A21;
fusion::at_c<1>(fusion::at_c<1>(A)) = A22;
fusion::at_c<0>(A)[_0] = A11;
fusion::at_c<0>(A)[_1] = A12;
fusion::at_c<1>(A)[_0] = A21;
fusion::at_c<1>(A)[_1] = A22;
printmatrix(std::cout,A11,"matrix A11","row",9,1);
printmatrix(std::cout,A12,"matrix A12","row",9,1);
......@@ -592,8 +596,8 @@ void test_MultiTypeBlockVector_MultiTypeBlockMatrix() {
loop.apply(x,b,r);
printvector(std::cout,fusion::at_c<0>(x),"solution x1","entry",11,9,1);
printvector(std::cout,fusion::at_c<1>(x),"solution x2","entry",11,9,1);
printvector(std::cout,x[_0],"solution x1","entry",11,9,1);
printvector(std::cout,x[_1],"solution x2","entry",11,9,1);
}
#endif
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment