Skip to content
Snippets Groups Projects
Commit b4961bb3 authored by Christian Engwer's avatar Christian Engwer
Browse files

* add fassigntest

* fix fs#652

[[Imported from SVN: r5724]]
parent abbcbb38
No related branches found
No related tags found
No related merge requests found
...@@ -28,12 +28,12 @@ namespace Dune { ...@@ -28,12 +28,12 @@ namespace Dune {
* overload operator <<= for FieldVector assignment from Dune::Zero * overload operator <<= for FieldVector assignment from Dune::Zero
*/ */
struct Zero { struct Zero {
Zero (int) {}; explicit Zero (int) {};
/** \brief Conversion operator to double */ /** \brief Conversion operator to double */
operator double () { return 0.0; } operator double () { return 0.0; }
/** \brief Conversion operator to int */ /** \brief Conversion operator to int */
operator int () { return 0; } operator int () { return 0; }
} zero = 0; } zero(0);
/** /**
* @brief Marker class for next row * @brief Marker class for next row
...@@ -41,8 +41,8 @@ namespace Dune { ...@@ -41,8 +41,8 @@ namespace Dune {
* overload operator <<= for FiledMatrix assignment * overload operator <<= for FiledMatrix assignment
*/ */
struct NextRow { struct NextRow {
NextRow (int) {}; explicit NextRow (int) {};
} nextRow = 0; } nextRow(0);
} // end empty namespace } // end empty namespace
...@@ -73,6 +73,7 @@ namespace Dune { ...@@ -73,6 +73,7 @@ namespace Dune {
FieldVector<T,s> & v; FieldVector<T,s> & v;
int c; int c;
bool temporary; bool temporary;
fvector_assigner();
public: public:
/*! @brief Copy Constructor */ /*! @brief Copy Constructor */
fvector_assigner(fvector_assigner & a) : v(a.v), c(a.c), temporary(false) fvector_assigner(fvector_assigner & a) : v(a.v), c(a.c), temporary(false)
...@@ -130,8 +131,8 @@ namespace Dune { ...@@ -130,8 +131,8 @@ namespace Dune {
* overload operator <<= for fvector assignment * overload operator <<= for fvector assignment
* from comma seperated list of values * from comma seperated list of values
*/ */
template <class T, int s> template <class T, class K, int s>
fvector_assigner<T,s> operator <<= (FieldVector<T,s> & v, const T & t) fvector_assigner<T,s> operator <<= (FieldVector<T,s> & v, const K & t)
{ {
return fvector_assigner<T,s>(v,true).append(t); return fvector_assigner<T,s>(v,true).append(t);
} }
...@@ -256,8 +257,8 @@ namespace Dune { ...@@ -256,8 +257,8 @@ namespace Dune {
* overload operator <<= for FieldMatrix assignment * overload operator <<= for FieldMatrix assignment
* from comma seperated list of values * from comma seperated list of values
*/ */
template <class T, int n, int m> template <class T, class K, int n, int m>
fmatrix_assigner<T,n,m> operator <<= (FieldMatrix<T,n,m> & v, const T & t) fmatrix_assigner<T,n,m> operator <<= (FieldMatrix<T,n,m> & v, const K & t)
{ {
return fmatrix_assigner<T,n,m>(v,true).append(t); return fmatrix_assigner<T,n,m>(v,true).append(t);
} }
......
Makefile
Makefile.in
.deps
.libs
semantic.cache
fvectortest
mpicollcomm
lrutest
parsetest
test-stack
arraylisttest
smartpointertest
float_cmp
bitsetvectortest
iteratorfacadetest
sllisttest
tuplestest
settest
fmatrixtest
poolallocatortest
*.gcda
*.gcno
gmon.out
gcdlcmtest
streamtest
exprtmpl
timing_xpr
timing_old
timing_flt
bigunsignedinttest
mpihelpertest
singletontest
utilitytest
testfassign_fail1
testfassign_fail2
testfassign_fail3
testfassign_fail4
testfassign_fail5
testfassign_fail6
testfassign1
testfassign2
testfassign3
testfassign4
smallobject
conversiontest
nullptr-test
blockbitfieldtest
deprtuplestest
fassigntest
...@@ -12,7 +12,7 @@ TESTPROGS = parsetest test-stack arraylisttest smartpointertest \ ...@@ -12,7 +12,7 @@ TESTPROGS = parsetest test-stack arraylisttest smartpointertest \
testfassign_fail1 testfassign_fail2 testfassign_fail3 \ testfassign_fail1 testfassign_fail2 testfassign_fail3 \
testfassign_fail4 testfassign_fail5 testfassign_fail6 \ testfassign_fail4 testfassign_fail5 testfassign_fail6 \
conversiontest bitsetvectortest deprtuplestest \ conversiontest bitsetvectortest deprtuplestest \
float_cmp float_cmp, fassigntest
# which tests to run # which tests to run
COMPILE_XFAIL=$(DUNE_COMMON_BIN)/xfail-compile-tests COMPILE_XFAIL=$(DUNE_COMMON_BIN)/xfail-compile-tests
...@@ -53,6 +53,8 @@ nullptr_test_fail_CPPFLAGS = -DFAIL ...@@ -53,6 +53,8 @@ nullptr_test_fail_CPPFLAGS = -DFAIL
static_assert_test_SOURCES = static_assert_test.cc static_assert_test_SOURCES = static_assert_test.cc
fassigntest_SOURCES = fassigntest.cc
bigunsignedinttest_SOURCES=bigunsignedinttest.cc bigunsignedinttest_SOURCES=bigunsignedinttest.cc
parsetest_SOURCES = parsetest.cc parsetest_SOURCES = parsetest.cc
......
// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
// vi: set et ts=4 sw=2 sts=2:
#include <config.h>
#include <iostream>
#include <dune/common/fvector.hh>
#include <dune/common/fassign.hh>
using namespace Dune;
int main(int argc, char** argv) try
{
Dune::FieldVector<double,3> pos;
pos <<= 1, 0, 0;
}
catch (Exception e) {
std::cout << e << std::endl;
}
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