12#include <deal.II/base/parameter_acceptor.h>
57 : ParameterAcceptor(subsection)
59 filename_ =
"ryujin.tif";
60 this->add_parameter(
"filename", filename_,
"GeoTIFF: image file to load");
62 transformation_ = {0., 0.01, 0., 0., 0., 0.01};
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 "
73 transformation_allow_out_of_bounds_queries_ =
false;
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.");
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.");
89 transformation_use_geotiff_origin_ =
false;
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.");
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");
104 height_scaling_ = 1.0;
105 this->add_parameter(
"height scaling",
107 "GeoTIFF: choose base point for height normalization "
108 "that is set to 0.: none, minimum, average, maximum");
110 const auto set_up = [&] {
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};
125 this->parse_parameters_call_back.connect(set_up);
129 DEAL_II_ALWAYS_INLINE
inline double
137 const double x = point[0];
139 if constexpr (dim >= 2)
141 const auto &[di, dj] = apply_inverse_transformation(x, y);
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;
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 <<
")"
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"));
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]);
180 std::cout <<
"index bounding box: (" << i_left <<
"," << j_left
181 <<
") and (" << i_right <<
"," << j_right <<
")" << std::endl;
185 const double i_ratio = std::fmod(di, 1.);
186 const double j_ratio = std::fmod(dj, 1.);
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]];
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]];
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;
197 return height_scaling_ * (v_jl * (1. - j_ratio) + v_jr * j_ratio);
211 void read_in_raster()
const
214 auto dataset_handle = GDALOpen(filename_.c_str(), GA_ReadOnly);
215 AssertThrow(dataset_handle,
216 dealii::ExcMessage(
"GDAL error: file not found"));
218 auto dataset = GDALDataset::FromHandle(dataset_handle);
219 Assert(dataset, dealii::ExcInternalError());
221 const auto driver = dataset->GetDriver();
223 driver_name_ = driver->GetMetadataItem(GDAL_DMD_LONGNAME);
224 if (dataset->GetProjectionRef() !=
nullptr)
225 driver_projection_ = dataset->GetProjectionRef();
230 dataset->GetRasterCount() == 1,
232 "GDAL driver error: currently we only support one raster"));
234 const auto raster_band = dataset->GetRasterBand(1);
236 AssertThrow(dataset->GetRasterXSize() == raster_band->GetXSize() &&
237 dataset->GetRasterYSize() == raster_band->GetYSize(),
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."));
251 raster_offset_ = {0, 0};
252 raster_size_ = {dataset->GetRasterXSize(), dataset->GetRasterYSize()};
254 raster_.resize(raster_size_[0] * raster_size_[1]);
255 const auto error_code = raster_band->RasterIO(
268 AssertThrow(error_code == 0,
270 "GDAL driver error: error reading in geotiff file"));
282 if (transformation_use_geotiff_) {
284 dataset->GetGeoTransform(affine_transformation_.data()) == CE_None;
286 dealii::ExcMessage(
"GDAL driver error: no geo transform "
287 "present in geotiff file"));
289 affine_transformation_ = transformation_;
291 affine_transformation_[2] *= -1.;
292 affine_transformation_[5] *= -1.;
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];
318 inverse_affine_transformation_[0] = affine_transformation_[0];
319 inverse_affine_transformation_[3] = affine_transformation_[3];
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];
330 GDALClose(dataset_handle);
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;
355 shift = *std::min_element(std::begin(raster_), std::end(raster_));
357 shift = *std::max_element(std::begin(raster_), std::end(raster_));
360 dealii::ExcInternalError());
361 const auto sum = std::reduce(std::begin(raster_), std::end(raster_));
362 shift = sum / raster_.size();
365 std::for_each(std::begin(raster_),
367 [&](
auto &element) { element -= shift; });
371 static constexpr auto message =
372 "ryujin has to be configured with GDAL support in order to read in "
374 AssertThrow(
false, dealii::ExcMessage(message));
380 DEAL_II_ALWAYS_INLINE
inline std::array<double, 2>
381 apply_transformation(
const double i,
const double j)
const
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;
390 DEAL_II_ALWAYS_INLINE
inline std::array<double, 2>
391 apply_inverse_transformation(
const double x,
const double y)
const
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]);
401 std::string filename_;
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_;
408 double height_scaling_;
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_;
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)