--- /dev/null
+Fixed: FullMatrix::print_formatted() and LAPACKFullMatrix::print_formatted
+would previously print NaN values as zero. This is now fixed.
+<br>
+(Simon Sticko, Martin Kronbichler, 2017/06/08)
for (size_type i=0; i<m(); ++i)
{
for (size_type j=0; j<n(); ++j)
- if (std::abs((*this)(i,j)) > threshold)
+ // we might have complex numbers, so use abs also to check for nan
+ // since there is no isnan on complex numbers
+ if (std::isnan(std::abs((*this)(i,j))))
+ out << std::setw(width) << (*this)(i,j) << ' ';
+ else if (std::abs((*this)(i,j)) > threshold)
out << std::setw(width)
<< (*this)(i,j) * number(denominator) << ' ';
else
for (size_type i=0; i<this->n_rows(); ++i)
{
for (size_type j=0; j<this->n_cols(); ++j)
- if (std::fabs(this->el(i,j)) > threshold)
+ // we might have complex numbers, so use abs also to check for nan
+ // since there is no isnan on complex numbers
+ if (std::isnan(std::abs((*this)(i,j))))
+ out << std::setw(width) << (*this)(i,j) << ' ';
+ else if (std::fabs(this->el(i,j)) > threshold)
out << std::setw(width)
<< this->el(i,j) * denominator << ' ';
else
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2017 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+
+// Check FullMatrix::print_formatted on NaN entry
+
+#include "../tests.h"
+
+#include <deal.II/lac/full_matrix.h>
+
+#include <iomanip>
+#include <fstream>
+
+int main()
+{
+ initlog();
+
+ FullMatrix<double> matrix(2,2);
+ matrix(0,0) = std::numeric_limits<double>::quiet_NaN();
+ matrix(0,1) = std::numeric_limits<double>::infinity();
+ matrix(1,1) = -std::numeric_limits<double>::infinity();
+
+ deallog << "Using print" << std::endl;
+ matrix.print(deallog.get_file_stream());
+
+ deallog << "Using print_formatted" << std::endl;
+ matrix.print_formatted(deallog.get_file_stream(), 3, true, 0, "0");
+
+ return 0;
+}
--- /dev/null
+
+DEAL::Using print
+nan inf
+0.0 -inf
+DEAL::Using print_formatted
+nan inf
+0 -inf
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2017 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+
+// Check LAPACKFullMatrix::print_formatted on NaN entry
+
+#include "../tests.h"
+
+#include <deal.II/lac/lapack_full_matrix.h>
+
+#include <iomanip>
+#include <fstream>
+
+int main()
+{
+ initlog();
+
+ LAPACKFullMatrix<double> matrix(2,2);
+ matrix(0,0) = std::numeric_limits<double>::quiet_NaN();
+ matrix(0,1) = std::numeric_limits<double>::infinity();
+ matrix(1,1) = -std::numeric_limits<double>::infinity();
+
+ deallog << "Using print_formatted" << std::endl;
+ matrix.print_formatted(deallog.get_file_stream(), 3, true, 0, "0");
+
+ return 0;
+}
--- /dev/null
+
+DEAL::Using print_formatted
+nan inf
+0 -inf