diff --git a/CHANGELOG.md b/CHANGELOG.md index ccd9042f05c69fbf5b6cd525e3492064e2f7d454..76b43cbf04ba4b35e200c232a296bb2ec86813dd 100644 --- a/CHANGELOG.md +++ b/CHANGELOG.md @@ -1,3 +1,22 @@ +Differences Between DuMu<sup>x</sup> 3.7 and DuMu<sup>x</sup> 3.6 +============================================= + +### General changes / structure + + +### Improvements and Enhancements + +- __IO/RasterImageWriter__: A tool now exists for writing `.pbm` and `.pgm` image files. + +### Immediate interface changes not allowing/requiring a deprecation period: + + +### Deprecated properties/classes/functions/files, to be removed after 3.7: + + +### New experimental features (possibly subject to backwards-incompatible changes in the future) + + Differences Between DuMu<sup>x</sup> 3.6 and DuMu<sup>x</sup> 3.5 ============================================= diff --git a/dumux/io/grid/gridmanager_sub.hh b/dumux/io/grid/gridmanager_sub.hh index ae52f3ac7538474a051689c517531a8056cdf63d..2c47f037ac91c5b741d56b4af81bbcc07e9e1baf 100644 --- a/dumux/io/grid/gridmanager_sub.hh +++ b/dumux/io/grid/gridmanager_sub.hh @@ -34,6 +34,7 @@ #if HAVE_DUNE_SUBGRID #include <dune/subgrid/subgrid.hh> #include <dumux/io/rasterimagereader.hh> +#include <dumux/io/rasterimagewriter.hh> #endif #ifndef DUMUX_IO_GRID_MANAGER_HH diff --git a/dumux/io/rasterimagedata.hh b/dumux/io/rasterimagedata.hh new file mode 100644 index 0000000000000000000000000000000000000000..b8b9a0ffea12fde112156938ef3339d4e623d238 --- /dev/null +++ b/dumux/io/rasterimagedata.hh @@ -0,0 +1,96 @@ +// -*- 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 3 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 InputOutput + * \brief A data class for raster image information + */ +#ifndef DUMUX_RASTER_IMAGE_DATA_HH +#define DUMUX_RASTER_IMAGE_DATA_HH + +#include <string> +#include <vector> +#include <fstream> + +namespace Dumux::Detail::RasterImageData { + +/*! + * \brief A struct that holds all information of the image format. + */ +struct Format +{ + std::string magicNumber; + std::string type; + std::string encoding; +}; + +/*! + * \brief A struct that contains all header data of the image. + */ +struct HeaderData +{ + Format format; + std::size_t nCols; + std::size_t nRows; + std::size_t maxValue = 1; +}; + +/*! + * \brief The return type of the reading functions. + * Holds the actual pixel values and the header data. + */ +template<class T> +class Result : private std::vector<T> +{ + using Parent = std::vector<T>; +public: + Result() = delete; + + /*! + * \brief Construct from data and header by copy + */ + Result(const std::vector<T>& data, const HeaderData& header) + : Parent(data) + , header_(header) + {} + + /*! + * \brief Construct from data and header by move + */ + Result(std::vector<T>&& data, HeaderData&& header) + : Parent(std::move(data)) + , header_(std::move(header)) + {} + + //! Returns the header data. + const HeaderData& header() const { return header_; } + + // expose some methods of std::vector + using Parent::operator[]; + using Parent::begin; + using Parent::end; + using Parent::size; + +private: + HeaderData header_; +}; + +} // end namespace Dumux::Detail::RasterImageData + +#endif diff --git a/dumux/io/rasterimagereader.hh b/dumux/io/rasterimagereader.hh index 341a14bfa434de1bfc4a675c4fef94121209418b..038377f3bfb952651d4b092e693541f03eba0017 100644 --- a/dumux/io/rasterimagereader.hh +++ b/dumux/io/rasterimagereader.hh @@ -36,6 +36,7 @@ #include <dune/common/exceptions.hh> #include <dumux/common/stringutilities.hh> +#include <dumux/io/rasterimagedata.hh> namespace Dumux { @@ -46,67 +47,13 @@ namespace Dumux { */ class NetPBMReader { -public: - /*! - * \brief A struct that holds all information of the image format. - */ - struct Format - { - std::string magicNumber; - std::string type; - std::string encoding; - }; + template<typename T> + using Result = Detail::RasterImageData::Result<T>; - /*! - * \brief A struct that contains all header data of the image. - */ - struct HeaderData - { - Format format; - std::size_t nCols; - std::size_t nRows; - std::size_t maxValue = 1; - }; + using HeaderData = Detail::RasterImageData::HeaderData; + using Format = Detail::RasterImageData::Format; - /*! - * \brief The return type of the reading functions. - * Holds the actual pixel values and the header data. - */ - template<class T> - class Result : private std::vector<T> - { - using Parent = std::vector<T>; - public: - Result() = delete; - - /*! - * \brief Construct from data and header by copy - */ - Result(const std::vector<T>& data, const HeaderData& header) - : Parent(data) - , header_(header) - {} - - /*! - * \brief Construct from data and header by move - */ - Result(std::vector<T>&& data, HeaderData&& header) - : Parent(std::move(data)) - , header_(std::move(header)) - {} - - //! Returns the header data. - const HeaderData& header() const { return header_; } - - // expose some methods of std::vector - using Parent::operator[]; - using Parent::begin; - using Parent::end; - using Parent::size; - - private: - HeaderData header_; - }; +public: /*! * \brief A helper function to retrieve the format from tokens of the file's first line. diff --git a/dumux/io/rasterimagewriter.hh b/dumux/io/rasterimagewriter.hh new file mode 100644 index 0000000000000000000000000000000000000000..44fd9caea66454366559043881e913bdd572056b --- /dev/null +++ b/dumux/io/rasterimagewriter.hh @@ -0,0 +1,136 @@ +// -*- 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 3 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 InputOutput + * \brief A simple writer class for raster images. + */ +#ifndef DUMUX_RASTER_IMAGE_WRITER_HH +#define DUMUX_RASTER_IMAGE_WRITER_HH + +#include <cassert> +#include <string> +#include <vector> +#include <fstream> +#include <sstream> +#include <algorithm> +#include <map> +#include <iterator> +#include <iostream> + +#include <dune/common/exceptions.hh> +#include <dumux/common/stringutilities.hh> +#include <dumux/io/rasterimagedata.hh> + +namespace Dumux { + +/*! + * \ingroup InputOutput + * \brief A simple reader class for the Netpbm format (https://en.wikipedia.org/wiki/Netpbm_format). + * So far, only black and white (*.pbm) and grayscale (*pgm) images are supported. + */ +class NetPBMWriter +{ + template<typename T> + using Result = Detail::RasterImageData::Result<T>; + + using HeaderData = Detail::RasterImageData::HeaderData; + +public: + + template<class ValueType> + static void write(const std::string& writeFileName, + Result<ValueType>& img, + const bool useDuneGridOrdering = true) + { + // Pass the image to the writer + writeRasterImageFile_(writeFileName, img, useDuneGridOrdering); + } + + template<class ValueType> + static void write(const std::string& writeFileName, + const std::size_t& nCols, + const std::size_t& nRows, + const std::string& magicNumber, + const std::string& type, + const std::string& encoding, + const std::vector<ValueType>& img, + const bool useDuneGridOrdering = true) + { + // Fill header data and collect img data + HeaderData headerData; + headerData.nCols = nCols; + headerData.nRows = nRows; + headerData.format.magicNumber = magicNumber; + Result<ValueType> result(std::move(img), std::move(headerData)); + + writeRasterImageFile_(writeFileName, result, useDuneGridOrdering); + } + + + /*! + * \brief Change the ordering of the pixels according + * to Dune's convention, shifting the origin from upper left to lower left. + * + * \param result The image's pixel values ordered from top to bottom. + */ + template<class T> + static void applyDuneGridOrdering(Result<T>& result) + { + auto tmp = result; + for (std::size_t i = 0; i < result.size(); i += result.header().nCols) + std::swap_ranges((result.begin() + i), (result.begin() + i + result.header().nCols), (tmp.end() - i - result.header().nCols)); + } + +private: + + template <class T> + static void writeRasterImageFile_(const std::string& writeFileName, + Result<T>& result, + const bool useDuneGridOrdering = true) + { + // This will reverse any reordering + if (useDuneGridOrdering) + applyDuneGridOrdering(result); + + // Write the corrected image and header information to a file + std::ofstream outfile(writeFileName, std::ios::trunc); + outfile << result.header().format.magicNumber << "\n"; + outfile << result.header().nCols << " " << result.header().nRows << "\n"; + if ((result.header().format.magicNumber == "P2") || (result.header().format.magicNumber == "P5")) + { + for (int i = 0; i < result.size(); i++) + outfile << result[i] << "\n"; + } + else + { + for (int i = 0; i < result.size(); i++) + { + if (i % result.header().nCols == 0) // wrap once per row + outfile << "\n"; + outfile << result[i]; + } + } + } + +}; + +} // namespace Dumux + +#endif diff --git a/test/io/rasterimagereader/test_rasterimagereader.cc b/test/io/rasterimagereader/test_rasterimagereader.cc index 3853f420e7e4e84b62dd0bc2de02e680973cb668..a125d4c69f4307b2707da3a6f5b39a54bb5bb5fb 100644 --- a/test/io/rasterimagereader/test_rasterimagereader.cc +++ b/test/io/rasterimagereader/test_rasterimagereader.cc @@ -26,6 +26,7 @@ #include <iostream> #include <dumux/io/rasterimagereader.hh> +#include <dumux/io/rasterimagewriter.hh> #include <dumux/io/container.hh> //////////////////////// @@ -41,6 +42,13 @@ int main(int argc, char** argv) // read an ASCII image without applying Dune's cell ordering (origin at upper left) auto blackAndWhiteImageASCII = NetPBMReader::readPBM("blackwhite_j.pbm", false); + auto blackAndWhiteImageASCIIOutput = blackAndWhiteImageASCII; + for (int i = 0; i < blackAndWhiteImageASCII.size(); i++) + { + if (i < 5) + blackAndWhiteImageASCIIOutput[i] = 1; + } + NetPBMWriter::write("blackwhite_j_new.pbm", blackAndWhiteImageASCIIOutput, false); // create a 2D array to show the image in the terminal std::vector<std::vector<bool>> printableBlackAndWhiteImage; @@ -87,7 +95,6 @@ int main(int argc, char** argv) std::cout << "Reading black/white (binary) failed" << std::endl; return 1; } - std::cout << std::endl; // test file where the dimensions are given in the same line as the magic number and comments are present @@ -125,6 +132,14 @@ int main(int argc, char** argv) // read an ASCII image without applying Dune's cell ordering (origin at upper left) auto grayScaleImageASCII = NetPBMReader::template readPGM<std::size_t>("grayscale_j.pgm", false); + auto grayScaleImageASCIIOutput = grayScaleImageASCII; + for (int i = 0; i < grayScaleImageASCII.size(); i++) + { + if (i < 5) + grayScaleImageASCIIOutput[i] = 125; + } + NetPBMWriter::write("grayscale_j_new.pgm", grayScaleImageASCIIOutput, false); + std::vector<std::vector<std::uint8_t>> printableGrayScaleImage; printableGrayScaleImage.resize(grayScaleImageASCII.header().nRows, std::vector<std::uint8_t>(grayScaleImageASCII.header().nCols)); NetPBMReader::fillImage(printableGrayScaleImage, grayScaleImageASCII);