Correctly compute the minimum and maximum values.

As there are no functors for computing the minimum and maximum,
we convert the std::max and std::min function pointers to
functors (which is not really nice.) Previously we were somehow
tricked into using std::greater and std::less, which of course do
return true or false and not what we need. Additionally, do more
excessive testing with different ranges.
This commit is contained in:
Markus Blatt 2015-01-23 20:48:53 +01:00
parent 6fe660a3d5
commit aac4cb7d66
2 changed files with 62 additions and 17 deletions

View File

@ -35,9 +35,10 @@
#include <dune/common/parallel/communicator.hh>
#include <dune/common/enumset.hh>
#include<algorithm>
#include<limits>
#include<type_traits>
#include <algorithm>
#include <functional>
#include <limits>
#include <type_traits>
namespace Opm
{
@ -164,6 +165,9 @@ public:
/// This function can either be used with a container, an operator, and an initial value
/// to compute a reduction. Or with tuples of them to compute multiple reductions with only
/// one global communication.
/// The possible functors needed can be constructed with Opm::Reduction::makeGlobalMaxFunctor(),
/// Opm::Reduction::makeGlobalMinFunctor(), and
/// Opm::Reduction::makeGlobalSumFunctor().
/// \tparam type of the container or the tuple of containers.
/// \tparam tyoe of the operator or a tuple of operators, examples are e.g.
/// Reduction::MaskIDOperator, Reduction::MaskToMinOperator,
@ -369,6 +373,9 @@ private:
template<typename BinaryOperator>
struct MaskToMinOperator
{
MaskToMinOperator(BinaryOperator b)
: b_(b)
{}
/// \brief Apply the underlying binary operator according to the mask.
///
/// If mask is 0 then t2 will be substituted by the lowest value,
@ -420,6 +427,9 @@ private:
template<typename BinaryOperator>
struct MaskToMaxOperator
{
MaskToMaxOperator(BinaryOperator b)
: b_(b)
{}
/// \brief Apply the underlying binary operator according to the mask.
///
/// If mask is 0 then t2 will be substituted by the maximum value,
@ -450,6 +460,37 @@ private:
private:
BinaryOperator b_;
};
/// \brief Create a functor for computing a global sum.
///
/// To be used with ParallelISTLInformation::computeReduction.
template<class T>
MaskIDOperator<std::plus<T> >
makeGlobalSumFunctor()
{
return MaskIDOperator<std::plus<T> >();
}
/// \brief Create a functor for computing a global maximum.
///
/// To be used with ParallelISTLInformation::computeReduction.
template<class T>
MaskToMinOperator<std::pointer_to_binary_function<const T&,const T&,const T&> >
makeGlobalMaxFunctor()
{
return MaskToMinOperator<std::pointer_to_binary_function<const T&,const T&,const T&> >
(std::pointer_to_binary_function<const T&,const T&,const T&>
((const T&(*)(const T&, const T&))std::max<T>));
}
/// \brief Create a functor for computing a global minimum.
///
/// To be used with ParallelISTLInformation::computeReduction.
template<class T>
MaskToMaxOperator<std::pointer_to_binary_function<const T&,const T&,const T&> >
makeGlobalMinFunctor()
{
return MaskToMaxOperator<std::pointer_to_binary_function<const T&,const T&,const T&> >
(std::pointer_to_binary_function<const T&,const T&,const T&>
((const T&(*)(const T&, const T&))std::min<T>));
}
} // end namespace Reduction
} // end namespace Opm

View File

@ -36,9 +36,10 @@
#include <opm/core/linalg/ParallelIstlInformation.hpp>
#include <functional>
#ifdef HAVE_DUNE_ISTL
BOOST_AUTO_TEST_CASE(tupleReductionTest)
void runSumMaxMinTest(const int offset)
{
int N=100;
const int N=100;
int start, end, istart, iend;
std::tie(start,istart,iend,end) = computeRegions(N);
Opm::ParallelISTLInformation comm(MPI_COMM_WORLD);
@ -46,17 +47,24 @@ BOOST_AUTO_TEST_CASE(tupleReductionTest)
std::vector<int> x(end-start);
assert(comm.indexSet()->size()==x.size());
for(auto i=comm.indexSet()->begin(), iend=comm.indexSet()->end(); i!=iend; ++i)
x[i->local()]=i->global();
x[i->local()]=i->global()+offset;
auto containers = std::make_tuple(x, x, x);
auto operators = std::make_tuple(Opm::Reduction::MaskIDOperator<std::plus<int> >(),
Opm::Reduction::MaskToMinOperator<std::greater<int> >(),
Opm::Reduction::MaskToMaxOperator<std::less< int> >());
auto operators = std::make_tuple(Opm::Reduction::makeGlobalSumFunctor<int>(),
Opm::Reduction::makeGlobalMaxFunctor<int>(),
Opm::Reduction::makeGlobalMinFunctor<int>());
auto values = std::make_tuple(0,0,100000);
auto oldvalues = values;
comm.computeReduction(containers,operators,values);
BOOST_CHECK(std::get<0>(values)==std::get<0>(oldvalues)+((N-1)*N)/2);
BOOST_CHECK(std::get<1>(values)==std::min(0, std::get<1>(oldvalues)));
BOOST_CHECK(std::get<2>(values)==std::max(N, std::get<2>(oldvalues)));
BOOST_CHECK(std::get<0>(values)==std::get<0>(oldvalues)+((N-1+2*offset)*N)/2);
BOOST_CHECK(std::get<1>(values)==std::max(N+offset-1, std::get<1>(oldvalues)));
BOOST_CHECK(std::get<2>(values)==std::min(offset, std::get<2>(oldvalues)));
}
BOOST_AUTO_TEST_CASE(tupleReductionTest)
{
runSumMaxMinTest(0);
runSumMaxMinTest(20);
runSumMaxMinTest(-20);
}
BOOST_AUTO_TEST_CASE(singleContainerReductionTest)
{
@ -69,13 +77,9 @@ BOOST_AUTO_TEST_CASE(singleContainerReductionTest)
assert(comm.indexSet()->size()==x.size());
for(auto i=comm.indexSet()->begin(), iend=comm.indexSet()->end(); i!=iend; ++i)
x[i->local()]=i->global();
auto containers = std::make_tuple(x, x, x);
auto operators = std::make_tuple(Opm::Reduction::MaskIDOperator<std::plus<int> >(),
Opm::Reduction::MaskToMinOperator<std::greater<int> >(),
Opm::Reduction::MaskToMaxOperator<std::less< int> >());
int value = 1;
int oldvalue = value;
comm.computeReduction(x,Opm::Reduction::MaskIDOperator<std::plus<int> >(),value);
comm.computeReduction(x,Opm::Reduction::makeGlobalSumFunctor<int>(),value);
BOOST_CHECK(value==oldvalue+((N-1)*N)/2);
}
#endif