ryujin 2.1.1 revision 0348cbb53a3e4b1da2a4c037e81f88f2d21ce219
geotiff_reader.h
Go to the documentation of this file.
1//
2// SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception
3// Copyright (C) 2024 - 2024 by the ryujin authors
4//
5
6#pragma once
7
9#include "lazy.h"
10#include "patterns_conversion.h"
11
12#include <deal.II/base/parameter_acceptor.h>
13#include <iomanip>
14
15#ifdef WITH_GDAL
16#include <cpl_conv.h>
17#include <gdal.h>
18#include <gdal_priv.h>
19#endif
20
21namespace ryujin
22{
25 none,
27 minimum,
29 average,
32 };
33} // namespace ryujin
34
35#ifndef DOXYGEN
41#endif
42
43namespace ryujin
44{
52 template <int dim>
53 class GeoTIFFReader : dealii::ParameterAcceptor
54 {
55 public:
56 GeoTIFFReader(const std::string subsection)
57 : ParameterAcceptor(subsection)
58 {
59 filename_ = "ryujin.tif";
60 this->add_parameter("filename", filename_, "GeoTIFF: image file to load");
61
62 transformation_ = {0., 0.01, 0., 0., 0., 0.01};
63 this->add_parameter(
64 "transformation",
65 transformation_,
66 "Array \"t[]\" describing an affine transformation between image "
67 "space (indices i and j from bottom left) and real coordinates (x "
68 "and y): x = t[0] + t[1] * i + t[2] * j, and y = t[3] + t[4] * i + "
69 "t[5] * j. (This transformation sets the origin of the image space "
70 "into the bottom left corner with index i to the right and index j "
71 "up)");
72
73 transformation_allow_out_of_bounds_queries_ = false;
74 this->add_parameter(
75 "transformation allow out of bounds queries",
76 transformation_allow_out_of_bounds_queries_,
77 "GeoTIFF: allow out-of-bounds queries. When set to true, the reader "
78 "returns constant extended values for coordinates that are outside "
79 "of the image range.");
80
81 transformation_use_geotiff_ = true;
82 this->add_parameter("transformation use geotiff",
83 transformation_use_geotiff_,
84 "GeoTIFF: read in transformation from GeoTIFF for "
85 "constructing the affine transformation. If set to "
86 "false the manually specified transformation "
87 "parameters will be used instead.");
88
89 transformation_use_geotiff_origin_ = false;
90 this->add_parameter(
91 "transformation use geotiff origin",
92 transformation_use_geotiff_origin_,
93 "GeoTIFF: read in affine shift (i.e., position of "
94 "lower left corner) from GeoTIFF for constructing "
95 "the affine transformation. If set to false the origin specified "
96 "in the transformation parameter will be used instead.");
97
98 height_normalization_ = HeightNormalization::minimum;
99 this->add_parameter("height normalization",
100 height_normalization_,
101 "GeoTIFF: choose base point for height normalization "
102 "that is set to 0.: none, minimum, average, maximum");
103
104 height_scaling_ = 1.0;
105 this->add_parameter("height scaling",
106 height_scaling_,
107 "GeoTIFF: choose base point for height normalization "
108 "that is set to 0.: none, minimum, average, maximum");
109
110 const auto set_up = [&] {
111#ifdef WITH_GDAL
112 /* Initial GDAL and reset all data: */
113 GDALAllRegister();
114 driver_name_ = "";
115 driver_projection_ = "";
116 affine_transformation_ = {0, 0, 0, 0, 0, 0};
117 inverse_affine_transformation_ = {0, 0, 0, 0, 0, 0};
118 raster_offset_ = {0, 0};
119 raster_size_ = {0, 0};
120 raster_.clear();
121#endif
122 };
123
124 set_up();
125 this->parse_parameters_call_back.connect(set_up);
126 }
127
128
129 DEAL_II_ALWAYS_INLINE inline double
130 compute_height(const dealii::Point<dim> &point) const
131 {
132 geotiff_guard_.ensure_initialized([&]() {
133 read_in_raster();
134 return true;
135 });
136
137 const double x = point[0];
138 double y = 0;
139 if constexpr (dim >= 2)
140 y = point[1];
141 const auto &[di, dj] = apply_inverse_transformation(x, y);
142
143 /* Check that we are in bounds: */
144
145 const bool in_bounds =
146 di > -1.0 && di < static_cast<double>(raster_size_[0]) + 1.0 &&
147 dj > -1.0 && dj < static_cast<double>(raster_size_[1]) + 1.0;
148
149#ifdef DEBUG_OUTPUT
150 if (!in_bounds) {
151 std::cout << std::setprecision(16);
152 std::cout << "Queried point out of bounds." << std::endl;
153 std::cout << "Point: " << point << std::endl;
154 std::cout << "Transformed coordinates: (" << di << "," << dj << ")"
155 << std::endl;
156 }
157#endif
158
159 AssertThrow(
160 transformation_allow_out_of_bounds_queries_ || in_bounds,
161 dealii::ExcMessage("Raster error: The requested point is outside "
162 "the image boundary of the geotiff file"));
163
164 /*
165 * Use a simple bilinear interpolation and ensure we never go below
166 * the minimum or above the maximum index.
167 */
168
169 const auto i_left = std::min(
170 std::max(static_cast<int>(std::floor(di)), 0), raster_size_[0]);
171 const auto i_right = std::min(
172 std::max(static_cast<int>(std::ceil(di)), 0), raster_size_[0]);
173 const auto j_left = std::min(
174 std::max(static_cast<int>(std::floor(dj)), 0), raster_size_[1]);
175 const auto j_right = std::min(
176 std::max(static_cast<int>(std::ceil(dj)), 0), raster_size_[1]);
177
178#ifdef DEBUG_OUTPUT
179 if (!in_bounds) {
180 std::cout << "index bounding box: (" << i_left << "," << j_left
181 << ") and (" << i_right << "," << j_right << ")" << std::endl;
182 }
183#endif
184
185 const double i_ratio = std::fmod(di, 1.);
186 const double j_ratio = std::fmod(dj, 1.);
187
188 const auto v_iljl = raster_[i_left + j_left * raster_size_[0]];
189 const auto v_irjl = raster_[i_right + j_left * raster_size_[0]];
190
191 const auto v_iljr = raster_[i_left + j_right * raster_size_[0]];
192 const auto v_irjr = raster_[i_right + j_right * raster_size_[0]];
193
194 const auto v_jl = v_iljl * (1. - i_ratio) + v_irjl * i_ratio;
195 const auto v_jr = v_iljr * (1. - i_ratio) + v_irjr * i_ratio;
196
197 return height_scaling_ * (v_jl * (1. - j_ratio) + v_jr * j_ratio);
198 }
199
200 /*
201 * Return the affine transformation information that is stored in the
202 * GeoTIFF image.
203 */
205
209
210 private:
211 void read_in_raster() const
212 {
213#ifdef WITH_GDAL
214 auto dataset_handle = GDALOpen(filename_.c_str(), GA_ReadOnly);
215 AssertThrow(dataset_handle,
216 dealii::ExcMessage("GDAL error: file not found"));
217
218 auto dataset = GDALDataset::FromHandle(dataset_handle);
219 Assert(dataset, dealii::ExcInternalError());
220
221 const auto driver = dataset->GetDriver();
222
223 driver_name_ = driver->GetMetadataItem(GDAL_DMD_LONGNAME);
224 if (dataset->GetProjectionRef() != nullptr)
225 driver_projection_ = dataset->GetProjectionRef();
226
227 /* For now we support one raster in the dataset: */
228
229 AssertThrow(
230 dataset->GetRasterCount() == 1,
231 dealii::ExcMessage(
232 "GDAL driver error: currently we only support one raster"));
233
234 const auto raster_band = dataset->GetRasterBand(1);
235
236 AssertThrow(dataset->GetRasterXSize() == raster_band->GetXSize() &&
237 dataset->GetRasterYSize() == raster_band->GetYSize(),
238 dealii::ExcMessage(
239 "GDAL driver error: the raster band has a different "
240 "dimension than the (global) raster dimension of the "
241 "geotiff image. This is not supported."));
242
243 /*
244 * FIXME: For now, we simply read in the entire geotiff file on
245 * each rank. In order to save memory for very large files it would
246 * be possible to create a bounding box for the all active cells of
247 * the triangulation and then only read in a small region for which
248 * we actually need data.
249 */
250
251 raster_offset_ = {0, 0};
252 raster_size_ = {dataset->GetRasterXSize(), dataset->GetRasterYSize()};
253
254 raster_.resize(raster_size_[0] * raster_size_[1]);
255 const auto error_code = raster_band->RasterIO(
256 GF_Read,
257 raster_offset_[0], /* x-offset of image region */
258 raster_offset_[1], /* y-offset of image region */
259 raster_size_[0], /* x-size of image region */
260 raster_size_[1], /* y-size of image region */
261 raster_.data(),
262 raster_size_[0], /* x-size of target buffer */
263 raster_size_[1], /* y-size of target buffer */
264 GDT_Float32,
265 0,
266 0);
267
268 AssertThrow(error_code == 0,
269 dealii::ExcMessage(
270 "GDAL driver error: error reading in geotiff file"));
271
272 /*
273 * Read in the affine transformation from the geotiff image.
274 *
275 * Note that this transformation differs from the one we use in the
276 * parameter file: GDAL uses typical image orientation: the origin
277 * of the dataset is in the "top left" corner (instead of bottom
278 * left) and the first (column) index goes to the right and the
279 * second (row) index goes down.
280 */
281
282 if (transformation_use_geotiff_) {
283 const auto success =
284 dataset->GetGeoTransform(affine_transformation_.data()) == CE_None;
285 AssertThrow(success,
286 dealii::ExcMessage("GDAL driver error: no geo transform "
287 "present in geotiff file"));
288 } else {
289 affine_transformation_ = transformation_;
290 /* Flip sign for j index (y-coordinate): */
291 affine_transformation_[2] *= -1.;
292 affine_transformation_[5] *= -1.;
293 }
294
295 /*
296 * Ensure that (i=0, j=raster_size[1]-1) corresponds to the user
297 * supplied (transformation_[0], transformation_[3]).
298 */
299 if (transformation_use_geotiff_ == false ||
300 transformation_use_geotiff_origin_ == false) {
301 const auto j_max = raster_size_[1] - 1;
302 affine_transformation_[0] =
303 transformation_[0] - j_max * affine_transformation_[2];
304 affine_transformation_[3] =
305 transformation_[3] - j_max * affine_transformation_[5];
306 }
307
308 /*
309 * Compute inverse transformation of
310 *
311 * x = t[0] + t[1] * i + t[2] * j, y = t[3] + t[4] * i + t[5] * j.
312 *
313 * namely:
314 *
315 * i = it[1] * (x - it[0]) + it[2] * (y - it[3])
316 * j = it[4] * (x - it[0]) + it[5] * (y - it[3])
317 */
318 inverse_affine_transformation_[0] = affine_transformation_[0];
319 inverse_affine_transformation_[3] = affine_transformation_[3];
320
321 const auto determinant =
322 affine_transformation_[1] * affine_transformation_[5] -
323 affine_transformation_[2] * affine_transformation_[4];
324 const auto inv = 1. / determinant;
325 inverse_affine_transformation_[1] = inv * affine_transformation_[5];
326 inverse_affine_transformation_[2] = inv * (-affine_transformation_[2]);
327 inverse_affine_transformation_[4] = inv * (-affine_transformation_[4]);
328 inverse_affine_transformation_[5] = inv * affine_transformation_[1];
329
330 GDALClose(dataset_handle);
331
332#ifdef DEBUG_OUTPUT
333 std::cout << std::setprecision(16);
334 std::cout << "GDAL: driver name = " << driver_name_;
335 std::cout << "\nGDAL: projection = " << driver_projection_;
336 std::cout << "\nGDAL: transformation =";
337 for (const auto &it : affine_transformation_)
338 std::cout << " " << it;
339 std::cout << "\nGDAL: inverse trafo =";
340 for (const auto &it : inverse_affine_transformation_)
341 std::cout << " " << it;
342 std::cout << "\nGDAL: raster offset =";
343 for (const auto &it : raster_offset_)
344 std::cout << " " << it;
345 std::cout << "\nGDAL: raster size =";
346 for (const auto &it : raster_size_)
347 std::cout << " " << it;
348 std::cout << std::endl;
349#endif
350
351 if (height_normalization_ != HeightNormalization::none) {
352 float shift = 0.;
353
354 if (height_normalization_ == HeightNormalization::minimum)
355 shift = *std::min_element(std::begin(raster_), std::end(raster_));
356 else if (height_normalization_ == HeightNormalization::maximum)
357 shift = *std::max_element(std::begin(raster_), std::end(raster_));
358 else {
359 Assert(height_normalization_ == HeightNormalization::average,
360 dealii::ExcInternalError());
361 const auto sum = std::reduce(std::begin(raster_), std::end(raster_));
362 shift = sum / raster_.size();
363 }
364
365 std::for_each(std::begin(raster_),
366 std::end(raster_),
367 [&](auto &element) { element -= shift; });
368 }
369
370#else
371 static constexpr auto message =
372 "ryujin has to be configured with GDAL support in order to read in "
373 "GeoTIFF images";
374 AssertThrow(false, dealii::ExcMessage(message));
375 __builtin_trap();
376#endif
377 }
378
379
380 DEAL_II_ALWAYS_INLINE inline std::array<double, 2>
381 apply_transformation(const double i, const double j) const
382 {
383 const auto &at = affine_transformation_;
384 const double x = at[0] + at[1] * i + at[2] * j;
385 const double y = at[3] + at[4] * i + at[5] * j;
386 return {x, y};
387 }
388
389
390 DEAL_II_ALWAYS_INLINE inline std::array<double, 2>
391 apply_inverse_transformation(const double x, const double y) const
392 {
393 const auto &iat = inverse_affine_transformation_;
394 const double i = iat[1] * (x - iat[0]) + iat[2] * (y - iat[3]);
395 const double j = iat[4] * (x - iat[0]) + iat[5] * (y - iat[3]);
396 return {i, j};
397 }
398
399 /* Runtime parameters: */
400
401 std::string filename_;
402
403 std::array<double, 6> transformation_;
404 bool transformation_allow_out_of_bounds_queries_;
405 bool transformation_use_geotiff_;
406 bool transformation_use_geotiff_origin_;
407 HeightNormalization height_normalization_;
408 double height_scaling_;
409
410 /* GDAL data structures: */
411
412 //
413 // We use a Lazy<t> wrapper for lazy initialization with efficient
414 // Schmidt's double checking. We simply ignore the bool type here.
415 //
416 Lazy<bool> geotiff_guard_;
417 mutable std::string driver_name_;
418 mutable std::string driver_projection_;
419 mutable std::array<double, 6> affine_transformation_;
420 mutable std::array<double, 6> inverse_affine_transformation_;
421 mutable std::array<int, 2> raster_offset_;
422 mutable std::array<int, 2> raster_size_;
423 mutable std::vector<float> raster_;
424 };
425} // namespace ryujin
auto & height_scaling() const
DEAL_II_ALWAYS_INLINE double compute_height(const dealii::Point< dim > &point) const
GeoTIFFReader(const std::string subsection)
auto & affine_transformation() const
auto & raster_offset() const
auto & raster_size() const
void ensure_initialized(const Callable &creator) const
#define ACCESSOR_READ_ONLY(member)
HeightNormalization