Skip to content
Snippets Groups Projects
Commit 1cc8a813 authored by Markus Blatt's avatar Markus Blatt
Browse files

Merged modification from trunk. Branch is up to date now.

[[Imported from SVN: r6942]]
No related branches found
No related tags found
No related merge requests found
......@@ -505,10 +505,10 @@ namespace Dune
//! infinity norm (row sum norm, how to generalize for blocks?)
typename FieldTraits<value_type>::real_type infinity_norm () const
{
ConstIterator it = begin();
if (it == end()) // empty matrix
if (size() == 0)
return 0.0;
ConstIterator it = begin();
typename remove_const< typename FieldTraits<value_type>::real_type >::type max = it->one_norm();
for (it = it + 1; it != end(); ++it)
max = std::max(max, it->one_norm());
......@@ -519,13 +519,13 @@ namespace Dune
//! simplified infinity norm (uses Manhattan norm for complex values)
typename FieldTraits<value_type>::real_type infinity_norm_real () const
{
ConstIterator it = begin();
if (it == end()) // empty matrix
if (size() == 0)
return 0.0;
ConstIterator it = begin();
typename remove_const< typename FieldTraits<value_type>::real_type >::type max = it->one_norm_real();
for (it = it + 1; it != end(); ++it)
max = std::max(max, it->one_norm());
max = std::max(max, it->one_norm_real());
return max;
}
......
......@@ -531,19 +531,29 @@ namespace Dune {
//! infinity norm (maximum of absolute values of entries)
typename FieldTraits<value_type>::real_type infinity_norm () const
{
typename FieldTraits<value_type>::real_type result( 0 );
for (size_type i=0; i<size(); i++)
result = std::max(result, std::abs((*this)[i]));
return result;
if (size() == 0)
return 0.0;
ConstIterator it = begin();
typename FieldTraits<value_type>::real_type max = std::abs(*it);
for (it = it + 1; it != end(); ++it)
max = std::max(max, std::abs(*it));
return max;
}
//! simplified infinity norm (uses Manhattan norm for complex values)
typename FieldTraits<value_type>::real_type infinity_norm_real () const
{
typename FieldTraits<value_type>::real_type result( 0 );
for (size_type i=0; i<size(); i++)
result = std::max(result, fvmeta::absreal((*this)[i]));
return result;
if (size() == 0)
return 0.0;
ConstIterator it = begin();
typename FieldTraits<value_type>::real_type max = fvmeta::absreal(*it);
for (it = it + 1; it != end(); ++it)
max = std::max(max, fvmeta::absreal(*it));
return max;
}
//===== sizes
......
......@@ -38,6 +38,8 @@ namespace Dune {
public:
/** \brief storage for key lists
*/
typedef std::vector<std::string> KeyVector;
/** \brief Create new empty ParameterTree
......
......@@ -10,6 +10,8 @@
#include <iostream>
#include <algorithm>
#include <vector>
#include <cassert>
#include <complex>
using namespace Dune;
......@@ -496,9 +498,50 @@ void test_invert ()
A.invert();
}
// Make sure that a matrix with only NaN entries has norm NaN.
// Prior to r6819, the infinity_norm would be zero; see also FS #1147.
void
test_nan()
{
double mynan = 0.0/0.0;
Dune::FieldMatrix<double, 2, 2> m2(mynan);
assert(std::isnan(m2.infinity_norm()));
assert(std::isnan(m2.frobenius_norm()));
assert(std::isnan(m2.frobenius_norm2()));
Dune::FieldMatrix<double, 0, 2> m02(mynan);
assert(0.0 == m02.infinity_norm());
assert(0.0 == m02.frobenius_norm());
assert(0.0 == m02.frobenius_norm2());
Dune::FieldMatrix<double, 2, 0> m20(mynan);
assert(0.0 == m20.infinity_norm());
assert(0.0 == m20.frobenius_norm());
assert(0.0 == m20.frobenius_norm2());
}
// The computation of infinity_norm_real() was flawed from r6819 on
// until r6915.
void
test_infinity_norms()
{
std::complex<double> threefour(3.0, -4.0);
std::complex<double> eightsix(8.0, -6.0);
Dune::FieldMatrix<std::complex<double>, 2, 2> m;
m[0] = threefour;
m[1] = eightsix;
assert(std::abs(m.infinity_norm() -20.0) < 1e-10); // max(5+5, 10+10)
assert(std::abs(m.infinity_norm_real()-28.0) < 1e-10); // max(7+7, 14+14)
}
int main()
{
try {
test_nan();
test_infinity_norms();
// test 1 x 1 matrices
test_matrix<float, 1, 1>();
ScalarOperatorTest<float>();
......
......@@ -11,6 +11,7 @@
#include <iostream>
#include <complex>
#include <typeinfo>
#include <cassert>
using Dune::FieldVector;
......@@ -308,12 +309,48 @@ public:
}
};
// Make sure that a vector with only NaN entries has norm NaN.
// Prior to r6914, the infinity_norm would be zero; see also FS #1147.
void
test_nan()
{
double mynan = 0.0/0.0;
Dune::FieldVector<double, 2> v2(mynan);
assert(std::isnan(v2.infinity_norm()));
assert(std::isnan(v2.one_norm()));
assert(std::isnan(v2.two_norm()));
assert(std::isnan(v2.two_norm2()));
Dune::FieldVector<double, 0> v0(mynan);
assert(0.0 == v0.infinity_norm());
assert(0.0 == v0.one_norm());
assert(0.0 == v0.two_norm());
assert(0.0 == v0.two_norm2());
}
void
test_infinity_norms()
{
std::complex<double> threefour(3.0, -4.0);
std::complex<double> eightsix(8.0, -6.0);
Dune::FieldVector<std::complex<double>, 2> v;
v[0] = threefour;
v[1] = eightsix;
assert(std::abs(v.infinity_norm() -10.0) < 1e-10); // max(5,10)
assert(std::abs(v.infinity_norm_real()-14.0) < 1e-10); // max(7,14)
}
int main()
{
try {
FieldVectorTest<int, 3>();
FieldVectorTest<float, 3>();
FieldVectorTest<double, 3>();
test_nan();
test_infinity_norms();
} catch (Dune::Exception& e) {
std::cerr << e << std::endl;
return 1;
......
......@@ -68,6 +68,8 @@ AC_DEFUN([DUNE_COMMON_CHECKS],
AC_CHECK_LIB([m], [pow])
AC_CHECK_FUNCS([sqrt strchr])
AC_LANG_POP([C++])
AC_REQUIRE([DUNE_PATH_GMP])
])
AC_DEFUN([DUNE_COMMON_CHECK_MODULE],
......
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