version 3.8
rasterimagewriter.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_RASTER_IMAGE_WRITER_HH
13#define DUMUX_RASTER_IMAGE_WRITER_HH
14
15#include <cassert>
16#include <string>
17#include <vector>
18#include <fstream>
19#include <sstream>
20#include <algorithm>
21#include <map>
22#include <iterator>
23#include <iostream>
24
25#include <dune/common/exceptions.hh>
28
29namespace Dumux {
30
37{
38 template<typename T>
40
42
43public:
44
45 template<class ValueType>
46 static void write(const std::string& writeFileName,
48 const bool useDuneGridOrdering = true)
49 {
50 // Pass the image to the writer
51 writeRasterImageFile_(writeFileName, img, useDuneGridOrdering);
52 }
53
54 template<class ValueType>
55 static void write(const std::string& writeFileName,
56 const std::size_t& nCols,
57 const std::size_t& nRows,
58 const std::string& magicNumber,
59 const std::string& type,
60 const std::string& encoding,
61 const std::vector<ValueType>& img,
62 const bool useDuneGridOrdering = true)
63 {
64 // Fill header data and collect img data
65 HeaderData headerData;
66 headerData.nCols = nCols;
67 headerData.nRows = nRows;
68 headerData.format.magicNumber = magicNumber;
69 Result<ValueType> result(std::move(img), std::move(headerData));
70
71 writeRasterImageFile_(writeFileName, result, useDuneGridOrdering);
72 }
73
74
81 template<class T>
82 static void applyDuneGridOrdering(Result<T>& result)
83 {
84 auto tmp = result;
85 for (std::size_t i = 0; i < result.size(); i += result.header().nCols)
86 std::swap_ranges((result.begin() + i), (result.begin() + i + result.header().nCols), (tmp.end() - i - result.header().nCols));
87 }
88
89private:
90
91 template <class T>
92 static void writeRasterImageFile_(const std::string& writeFileName,
93 Result<T>& result,
94 const bool useDuneGridOrdering = true)
95 {
96 // This will reverse any reordering
97 if (useDuneGridOrdering)
99
100 // Write the corrected image and header information to a file
101 std::ofstream outfile(writeFileName, std::ios::trunc);
102 outfile << result.header().format.magicNumber << "\n";
103 outfile << result.header().nCols << " " << result.header().nRows << "\n";
104 if ((result.header().format.magicNumber == "P2") || (result.header().format.magicNumber == "P5"))
105 {
106 for (int i = 0; i < result.size(); i++)
107 outfile << result[i] << "\n";
108 }
109 else
110 {
111 for (int i = 0; i < result.size(); i++)
112 {
113 if (i % result.header().nCols == 0) // wrap once per row
114 outfile << "\n";
115 outfile << result[i];
116 }
117 }
118 }
119
120};
121
122} // namespace Dumux
123
124#endif
The return type of the reading functions. Holds the actual pixel values and the header data.
Definition: rasterimagedata.hh:48
const HeaderData & header() const
Returns the header data.
Definition: rasterimagedata.hh:70
A simple reader class for the Netpbm format (https://en.wikipedia.org/wiki/Netpbm_format)....
Definition: rasterimagewriter.hh:37
static void applyDuneGridOrdering(Result< T > &result)
Change the ordering of the pixels according to Dune's convention, shifting the origin from upper left...
Definition: rasterimagewriter.hh:82
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)
Definition: rasterimagewriter.hh:55
static void write(const std::string &writeFileName, Result< ValueType > &img, const bool useDuneGridOrdering=true)
Definition: rasterimagewriter.hh:46
Definition: adapt.hh:17
A data class for raster image information.
Helpers for working with strings.
std::string magicNumber
Definition: rasterimagedata.hh:26
A struct that contains all header data of the image.
Definition: rasterimagedata.hh:35
std::size_t nRows
Definition: rasterimagedata.hh:38
Format format
Definition: rasterimagedata.hh:36
std::size_t nCols
Definition: rasterimagedata.hh:37