Skip to content
Snippets Groups Projects
Commit 829065ed authored by Kilian Weishaupt's avatar Kilian Weishaupt
Browse files

Merge branch 'feature/partial' into 'master'

[common] Implement partial function for referencing parts of a tuple

See merge request !1353
parents 37d2845f 8b458f6a
No related branches found
No related tags found
1 merge request!1353[common] Implement partial function for referencing parts of a tuple
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set et ts=4 sw=4 sts=4:
/*****************************************************************************
* See the file COPYING for full copying permissions. *
* *
* This program is free software: you can redistribute it and/or modify *
* it under the terms of the GNU General Public License as published by *
* the Free Software Foundation, either version 2 of the License, or *
* (at your option) any later version. *
* *
* This program is distributed in the hope that it will be useful, *
* but WITHOUT ANY WARRANTY; without even the implied warranty of *
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
* GNU General Public License for more details. *
* *
* You should have received a copy of the GNU General Public License *
* along with this program. If not, see <http://www.gnu.org/licenses/>. *
*****************************************************************************/
/*!
* \file
* \ingroup Common
* \brief Get only parts of a container or tuple
*/
#ifndef DUMUX_COMMON_PARTIAL_HH
#define DUMUX_COMMON_PARTIAL_HH
#include <tuple>
#include <type_traits>
#include <dune/istl/multitypeblockvector.hh>
namespace Dumux {
/*!
* \brief a function to get a MultiTypeBlockVector with references to some entries of another MultiTypeBlockVector
* \param v a MultiTypeBlockVector
* \param indices the indices of the entries that should be referenced
*/
template<class ...Args, std::size_t ...i>
auto partial(Dune::MultiTypeBlockVector<Args...>& v, Dune::index_constant<i>... indices)
{
return Dune::MultiTypeBlockVector<std::add_lvalue_reference_t<std::decay_t<std::tuple_element_t<indices, std::tuple<Args...>>>>...>(v[indices]...);
}
/*!
* \brief a function to get a tuple with references to some entries of another tuple
* \param v a tuple
* \param indices a tuple of indices of the entries that should be referenced
*/
template<class ...Args, std::size_t ...i>
auto partial(std::tuple<Args...>& v, Dune::index_constant<i>... indices)
{
return std::tuple<std::add_lvalue_reference_t<std::decay_t<std::tuple_element_t<indices, std::tuple<Args...>>>>...>(std::get<indices>(v)...);
}
/*!
* \brief a function to get a MultiTypeBlockVector with references to some entries of another MultiTypeBlockVector
* \param t an std::tuple or Dune::MultiTypeBlockVector
* \param indices a tuple of indices of the entries that should be referenced
*/
template<class T, std::size_t ...i>
auto partial(T& t, std::tuple<Dune::index_constant<i>...> indices)
{
return partial(t, Dune::index_constant<i>{}...);
}
} // end namespace Dumux
#endif
......@@ -6,3 +6,5 @@ add_subdirectory(propertysystem)
add_subdirectory(spline)
add_subdirectory(timeloop)
add_subdirectory(typetraits)
dune_add_test(SOURCES test_partial.cc)
#include <config.h>
#include <iostream>
#include <tuple>
#include <dune/common/indices.hh>
#include <dune/common/classname.hh>
#include <dune/common/exceptions.hh>
#include <dune/common/float_cmp.hh>
#include <dune/common/fvector.hh>
#include <dune/istl/bvector.hh>
#include <dune/istl/multitypeblockvector.hh>
#include <dumux/common/partial.hh>
namespace Dumux {
template<template<class...> class T>
void runTest()
{
using Block1 = Dune::BlockVector<Dune::FieldVector<double, 1>>;
using Block2 = Dune::BlockVector<Dune::FieldVector<double, 3>>;
Block1 a, b;
Block2 c;
a = {1, 3, 4};
b = {4, 5, 6, 7};
c = {{1, 2, 3}, {1, 2, 3}};
using namespace Dune::Indices;
T<Block1, Block1, Block2> m;
std::cout << "Testing " << Dune::className(m) << '\n' << std::endl;
std::get<0>(m) = a;
std::get<1>(m) = b;
std::get<2>(m) = c;
auto p = partial(m, _0, _2);
if (!std::is_same<T<Block1&, Block2&>, std::decay_t<decltype(p)>>::value)
DUNE_THROW(Dune::Exception, "Dumux::partial() returned wrong type: " << Dune::className(p));
std::get<1>(p)[0][0] = 5.0;
if (!Dune::FloatCmp::eq(std::get<2>(m)[0][0], 5.0))
DUNE_THROW(Dune::Exception, "Modifying referenced partial vector failed! (m = " << std::get<2>(m)[0][0] << ", p = " << std::get<1>(p)[0][0] << ")");
}
} // end namespace Dumux
int main(int argc, char* argv[]) try
{
using namespace Dumux;
runTest<Dune::MultiTypeBlockVector>();
runTest<std::tuple>();
return 0;
}
catch (const Dune::Exception& e)
{
std::cout << e << std::endl;
return 1;
}
catch (...)
{
std::cout << "Unknown exception thrown!" << std::endl;
return 1;
}
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