Skip to content
Snippets Groups Projects
Commit 5304152a authored by Christoph Grüninger's avatar Christoph Grüninger
Browse files

[!735] Fix inconsistent chained DenseVector ctors/assignment.

Merge branch 'fix/fix-inconsistent-1d-dense-vector-conversions' into 'master'

ref:core/dune-common The commit also adds a test which fails without the other
changes in the commit.

Problem is that DenseVector and FieldVector\<\..., 1\> implement catch-all
assignment constructors and ctors. I.e., the 1d FieldVector provides CTOR
call-signatures and the DenseVector<C> provides and assignment operator

std::is_assignable<To, From\> std::is_convertible<From, To\>

This poses a problem when nesting vectors and trying to assign to them from
another DenseVector implementation.

The problem is that the implementation of the CTOR and the assignment
operation fails, but OTOH the FieldVector implementation uses
std::is_convertible. Now std::convertible will report "yes" but the body of
the assignment operation will fail nevertheless.

The commit fixes this by restricting the call-signature to the cases where the
implementation of the assignment and the CTOR will succeed but testing the
assignability of DV::value_type.

Supersedes [!562] as this is a rebased and squased version.

See merge request [!735]

  [!562]: gitlab.dune-project.org/NoneNone/merge_requests/562
  [!735]: gitlab.dune-project.org/core/dune-common/merge_requests/735
parents 0ce22beb bed4e242
No related branches found
No related tags found
1 merge request!735Fix inconsistent chained DenseVector ctors/assignment.
Pipeline #23073 passed
......@@ -279,7 +279,9 @@ namespace Dune {
public:
//! Assignment operator for other DenseVector of different type
template <typename W>
template <typename W,
std::enable_if_t<
std::is_assignable<value_type&, typename DenseVector<W>::value_type>::value, int> = 0>
derived_type& operator= (const DenseVector<W>& other)
{
assert(other.size() == size());
......
......@@ -318,14 +318,16 @@ namespace Dune {
template<typename T,
typename EnableIf = typename std::enable_if<
std::is_convertible<T, K>::value &&
! std::is_same<K, DenseVector<typename FieldTraits<T>::field_type>
>::value
! std::is_base_of<DenseVector<typename FieldTraits<T>::field_type>, K
>::value
>::type
>
FieldVector (const T& k) : _data(k) {}
//! Constructor from static vector of different type
template<class C>
template<class C,
std::enable_if_t<
std::is_assignable<K&, typename DenseVector<C>::value_type>::value, int> = 0>
FieldVector (const DenseVector<C> & x)
{
static_assert(((bool)IsFieldVectorSizeCorrect<C,1>::value), "FieldVectors do not match in dimension!");
......@@ -359,9 +361,9 @@ namespace Dune {
//! Assignment operator for scalar
template<typename T,
typename EnableIf = typename std::enable_if<
std::is_convertible<T, K>::value &&
! std::is_same<K, DenseVector<typename FieldTraits<T>::field_type>
>::value
std::is_assignable<K&, T>::value &&
! std::is_base_of<DenseVector<typename FieldTraits<T>::field_type>, K
>::value
>::type
>
inline FieldVector& operator= (const T& k)
......
......@@ -201,6 +201,10 @@ dune_add_test(SOURCES fvectortest.cc
LINK_LIBRARIES dunecommon
LABELS quick)
dune_add_test(SOURCES fvectorconversion1d.cc
LINK_LIBRARIES dunecommon
LABELS quick)
dune_add_test(SOURCES genericiterator_compile_fail.cc
EXPECT_COMPILE_FAIL
LABELS quick)
......
// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
// vi: set et ts=4 sw=2 sts=2:
#ifdef HAVE_CONFIG_H
#include "config.h"
#endif
#include <iostream>
#include <algorithm>
#include <dune/common/densevector.hh>
#include <dune/common/fvector.hh>
#include <dune/common/fmatrix.hh>
#include <dune/common/exceptions.hh>
#include <dune/common/unused.hh>
template<class Component, std::size_t Dim >
class MyVector;
namespace Dune
{
template<class Component, std::size_t Dim>
struct DenseMatVecTraits< MyVector<Component, Dim> >
{
using derived_type = MyVector<Component, Dim>;
using value_type = Component;
using size_type = std::size_t;
};
template<class Component, std::size_t Dim, int Size>
struct IsFieldVectorSizeCorrect<MyVector<Component, Dim>, Size>
: std::integral_constant<bool, Dim == Size>
{};
}
template<class Component, std::size_t Dim >
class MyVector : public Dune::DenseVector< MyVector<Component, Dim> >
{
public:
static constexpr std::size_t size () { return Dim; }
Component& operator[] ( std::size_t i ) { return data_; }
const Component& operator[] ( std::size_t i ) const { return data_; }
protected:
Component data_;
};
int main()
{
try
{
// Pure 1d case. Here OuterMV is assignable to MiddleFV as the the
// 1d FieldVector implements a type-case to the underlying
// field. This is expected behaviour.
{
using InnerFV = Dune::FieldVector<double, 1>;
using MiddleFV = Dune::FieldVector<InnerFV, 1>;
using OuterFV = Dune::FieldVector<MiddleFV, 1>;
using MiddleMV = MyVector<InnerFV, 1>;
using OuterMV = MyVector<MiddleMV, 1>;
MiddleFV mfv;
OuterMV mv;
OuterFV fv;
static_assert(std::is_convertible<OuterMV, OuterFV>::value,
"DenseVectors should be convertible.");
fv = mv;
static_assert(std::is_assignable<MiddleFV&, OuterMV>::value,
"Reduced assignability detected.");
mfv = mv;
}
// The following will trigger a problem in the DenseVector
// operator=() which can be cured by first checking whether the
// value_types are assignable.
{
using InnerFV = Dune::FieldVector<double, 2>;
using MiddleFV = Dune::FieldVector<InnerFV, 1>;
using OuterFV = Dune::FieldVector<MiddleFV, 1>;
using MiddleMV = MyVector<InnerFV, 1>;
using OuterMV = MyVector<MiddleMV, 1>;
MiddleFV mfv;
OuterMV mv;
OuterFV fv;
static_assert(std::is_convertible<OuterMV, OuterFV>::value,
"DenseVectors should be convertible.");
fv = mv;
static_assert(!std::is_assignable<MiddleFV&, OuterMV>::value,
"Inconsistent assignability detected.");
// mfv = mv; // <- this would not compile dispite the assignability check.
}
{
using InnerFV = Dune::FieldMatrix<double, 2, 2>;
using MiddleFV = Dune::FieldVector<InnerFV, 1>;
using OuterFV = Dune::FieldVector<MiddleFV, 1>;
using MiddleMV = MyVector<InnerFV, 1>;
using OuterMV = MyVector<MiddleMV, 1>;
MiddleFV mfv;
OuterMV mv;
OuterFV fv;
static_assert(std::is_assignable<OuterFV, OuterMV>::value,
"DenseVectors should be assignable.");
fv = mv;
static_assert(!std::is_assignable<MiddleFV&, OuterMV>::value,
"Inconsistent assignability detected.");
// mfv = mv; // <- this would not compile dispite the assignability check.
}
return 0;
} catch (Dune::Exception& e) {
std::cerr << e << std::endl;
return 1;
} catch (...) {
std::cerr << "Generic exception!" << std::endl;
return 2;
}
}
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