PRISMS-PF Manual
Loading...
Searching...
No Matches
read_binary.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: © 2025 PRISMS Center at the University of Michigan
2// SPDX-License-Identifier: GNU Lesser General Public Version 2.1
3
4#pragma once
5
6#include <deal.II/base/exceptions.h>
7#include <deal.II/base/point.h>
8#include <deal.II/lac/vector.h>
9
10#include <prismspf/core/types.h>
11
13
15
17
18#include <bit>
19#include <fstream>
20#include <iostream>
21#include <ranges>
22#include <string>
23#include <utility>
24
26
30template <unsigned int dim, typename number>
31class ReadBinary : public ReadFieldBase<dim, number>
32{
33public:
37 ReadBinary(const InitialConditionFile &_ic_file,
38 const SpatialDiscretization<dim> &_spatial_discretization);
39
43 void
44 print_file() override;
45
49
50 static void
51 write_file(const std::vector<number> &data, const InitialConditionFile &ic_file);
52
56 number
57 get_scalar_value(const dealii::Point<dim> &point,
58 const std::string &scalar_name) override;
59
63 dealii::Vector<number>
64 get_vector_value(const dealii::Point<dim> &point,
65 const std::string &vector_name) override;
66
67private:
72 void
74
78 dealii::Vector<number>
79 get_value(const dealii::types::global_dof_index index, const unsigned int n_components);
80
90 dealii::Vector<number>
91 interpolate(const dealii::Point<dim> &point, const unsigned int n_components);
92
99 dealii::types::global_dof_index n_points = 1;
100
104 dealii::types::global_dof_index n_values = 0;
105
109 std::vector<number> data;
110};
111
112template <unsigned int dim, typename number>
114 const InitialConditionFile &_ic_file,
115 const SpatialDiscretization<dim> &_spatial_discretization)
116 : ReadFieldBase<dim, number>(_ic_file, _spatial_discretization)
117{
118 // Make sure the dataset format is correct
119 AssertThrow(this->ic_file.dataset_format == DataFormatType::FlatBinary,
120 dealii::ExcMessage("Dataset format must be FlatBinary"));
121
122 // Check that only one field is being read in
123 AssertThrow(this->ic_file.file_variable_names.size() == 1 &&
124 this->ic_file.simulation_variable_names.size() == 1,
125 dealii::ExcMessage("Only one field can be read in from a binary file"));
126
127 // Make sure we have a rectangular domain
128 AssertThrow(_spatial_discretization.type == TriangulationType::Rectangular,
129 dealii::ExcMessage(
130 "Only rectangular domains are supported for binary input files"));
131
132 // Compute the total number of points in the binary file
133 for (unsigned int d : std::views::iota(0U, dim))
134 {
135 n_points *= this->ic_file.n_data_points[d];
136 }
137
138 // Check that the binary matches an expected size
140
141 // Read in the binary file
142 std::ifstream data_file(this->ic_file.filename, std::ios::binary);
143 AssertThrow(data_file,
144 dealii::ExcMessage("Could not open binary file: " +
145 this->ic_file.filename));
146
147 // Reserve space in the data vector
148 data.resize(n_values);
149
150 // Read in the data
151 data_file.read(reinterpret_cast<char *>(data.data()), n_values * sizeof(number));
152 data_file.close();
153}
154
155template <unsigned int dim, typename number>
156inline void
158{
159 // Grab the file size of the binary file in bytes
160 auto file_size = std::filesystem::file_size(this->ic_file.filename);
161
162 // Compute the expected size of the binary file. This is simply the number of points
163 // multiplied by the size of each point in bytes.
164 auto expected_size_scalar = static_cast<std::uintmax_t>(n_points * sizeof(number));
165 auto expected_size_vector = static_cast<std::uintmax_t>(dim * expected_size_scalar);
166
167 // Make sure expected size is not zero
168 AssertThrow(
169 expected_size_scalar != 0 && expected_size_vector != 0,
170 dealii::ExcMessage(
171 "Expected input array size is zero, check that the number of data points "
172 "in each used direction is set correctly in the input file for your binary file. "
173 "You likely have the number of data points set to zero in all directions."));
174 // Make sure the size matches for either a scalar or vector
175 AssertThrow(file_size == expected_size_scalar || file_size == expected_size_vector,
176 dealii::ExcMessage(
177 "Expected binary file size (" + std::to_string(expected_size_scalar) +
178 " bytes for scalar or " + std::to_string(expected_size_vector) +
179 " bytes for vector) does not match actual file size (" +
180 std::to_string(file_size) + " bytes)."));
181
182 // Set the number of values
183 n_values = file_size == expected_size_scalar ? n_points : dim * n_points;
184}
185
186template <unsigned int dim, typename number>
187inline void
188ReadBinary<dim, number>::write_file(const std::vector<number> &data,
190{
191 // Try to open the file
192 std::ofstream data_file(ic_file.filename, std::ios::binary);
193 AssertThrow(data_file,
194 dealii::ExcMessage("Could not open binary file: " + ic_file.filename));
195
196 // Write the data
197 data_file.write(reinterpret_cast<const char *>(data.data()),
198 data.size() * sizeof(number));
199 data_file.close();
200}
201
202template <unsigned int dim, typename number>
203inline dealii::Vector<number>
204ReadBinary<dim, number>::get_value(const dealii::types::global_dof_index index,
205 const unsigned int n_components)
206{
207 // Create vector to hold the value with the correct number of components
208 dealii::Vector<number> value(n_components);
209
210 // Check that the number of components matches what we expect from the number of values
211 Assert(n_values % n_components == 0,
212 dealii::ExcMessage("The number of components requested in the get_value call "
213 "does not match the number of values in the binary file."));
214
215 // Fill the value vector
216 for (unsigned int component : std::views::iota(0U, n_components))
217 {
218 Assert(((n_components * index) + component) < data.size(),
219 dealii::ExcMessage("Index out of bounds in ReadBinary::get_vals"));
220 value[component] = data[(n_components * index) + component];
221 }
222
223 return value;
224}
225
226template <unsigned int dim, typename number>
227inline dealii::Vector<number>
228ReadBinary<dim, number>::interpolate(const dealii::Point<dim> &point,
229 const unsigned int n_components)
230{
231 Assert(n_components == 1 || n_components == dim,
232 dealii::ExcMessage(
233 "Number of components for interpolation must be 1 (scalar) or dim (vector)"));
234
235 // Create a vector to hold the interpolated value
236 dealii::Vector<number> value(n_components);
237
238 // Get the spatial discretization
240
241 // Compute the spacing in each direction
242 std::array<number, dim> spacing;
243 for (unsigned int d : std::views::iota(0U, dim))
244 {
245 spacing[d] = spatial_discretization.rectangular_mesh.size[d] /
246 static_cast<number>(this->ic_file.n_data_points[d]);
247 }
248
249 // Compute the indices of the lower corner of the cell containing the point
250 std::array<dealii::types::global_dof_index, dim> lower_indices;
251 for (unsigned int d : std::views::iota(0U, dim))
252 {
253 lower_indices[d] =
254 static_cast<dealii::types::global_dof_index>(std::floor(point[d] / spacing[d]));
255 // Make sure we don't go out of bounds
256 if (lower_indices[d] >= this->ic_file.n_data_points[d] - 1)
257 {
258 lower_indices[d] = this->ic_file.n_data_points[d] - 2;
259 }
260 }
261
262 // Compute the weights for interpolation in each direction
263 std::array<number, dim> weights;
264 for (unsigned int d : std::views::iota(0U, dim))
265 {
266 weights[d] = (point[d] - lower_indices[d] * spacing[d]) / spacing[d];
267 }
268
269 // Perform multilinear interpolation based on the dimension
270 if constexpr (dim == 1)
271 {
272 // Here is the map of the nodes in 1D:
273 // 0 — 1
274
275 // Grab the index of the lower left corner of the cell
276 // cppcheck-suppress-begin containerOutOfBounds
277 auto lower_index = lower_indices[0];
278 // cppcheck-suppress-end containerOutOfBounds
279
280 // Get the indices of the two points in 1D
281 auto node_index_0 = lower_index;
282 auto node_index_1 = lower_index + 1;
283
284 // Get the values at these points
285 auto value_0 = get_value(node_index_0, n_components);
286 auto value_1 = get_value(node_index_1, n_components);
287
288 // Interpolate
289 for (unsigned int c : std::views::iota(0U, n_components))
290 {
291 // cppcheck-suppress-begin containerOutOfBounds
292 value[c] = (1.0 - weights[0]) * value_0[c] + weights[0] * value_1[c];
293 // cppcheck-suppress-end containerOutOfBounds
294 }
295 }
296 else if constexpr (dim == 2)
297 {
298 // Here is the map of the nodes in 2D:
299 // 01 — 11
300 // | |
301 // 00 — 10
302
303 // Grab the row length in the x-direction (0th direction)
304 auto row_length_0 = this->ic_file.n_data_points[0];
305
306 // Grab the index of the lower left corner of the cell
307 // cppcheck-suppress-begin containerOutOfBounds
308 auto lower_index = lower_indices[0] + (lower_indices[1] * row_length_0);
309 // cppcheck-suppress-end containerOutOfBounds
310
311 // Get the indices of the four points in 2D
312 auto node_index_00 = lower_index;
313 auto node_index_10 = lower_index + 1;
314
315 auto node_index_01 = lower_index + row_length_0;
316 auto node_index_11 = lower_index + row_length_0 + 1;
317
318 // Get the values at these points
319 auto value_00 = get_value(node_index_00, n_components);
320 auto value_10 = get_value(node_index_10, n_components);
321
322 auto value_01 = get_value(node_index_01, n_components);
323 auto value_11 = get_value(node_index_11, n_components);
324
325 // Interpolate
326 for (unsigned int c : std::views::iota(0U, n_components))
327 {
328 // cppcheck-suppress-begin containerOutOfBounds
329 value[c] = (1.0 - weights[0]) * (1.0 - weights[1]) * value_00[c] +
330 weights[0] * (1.0 - weights[1]) * value_10[c] +
331 (1.0 - weights[0]) * weights[1] * value_01[c] +
332 weights[0] * weights[1] * value_11[c];
333 // cppcheck-suppress-end containerOutOfBounds
334 }
335 }
336 else if constexpr (dim == 3)
337 {
338 // Here is the map of the nodes in 3D:
339 //
340 // 011 ———— 111
341 // / | / |
342 // / 001 / |
343 // 010 ——— 110 101
344 // | / | /
345 // | / | /
346 // 000 ———— 100
347
348 // Grab the row length in the x-direction (0th direction)
349 auto row_length_0 = this->ic_file.n_data_points[0];
350 // Grab the row length in the y-direction (1st direction)
351 auto row_length_1 = this->ic_file.n_data_points[1];
352
353 // Grab the index of the lower left corner of the cell
354 // cppcheck-suppress-begin containerOutOfBounds
355 auto lower_index = lower_indices[0] + (lower_indices[1] * row_length_0) +
356 (lower_indices[2] * row_length_0 * row_length_1);
357 // cppcheck-suppress-end containerOutOfBounds
358
359 // Get the indices of the eight points in 3D
360 auto node_index_000 = lower_index;
361 auto node_index_100 = lower_index + 1;
362
363 auto node_index_010 = lower_index + row_length_0;
364 auto node_index_110 = lower_index + row_length_0 + 1;
365
366 auto node_index_001 = lower_index + (row_length_0 * row_length_1);
367 auto node_index_101 = lower_index + (row_length_0 * row_length_1) + 1;
368
369 auto node_index_011 = lower_index + (row_length_0 * row_length_1) + row_length_0;
370 auto node_index_111 =
371 lower_index + (row_length_0 * row_length_1) + row_length_0 + 1;
372
373 // Get the values at these points
374 auto value_000 = get_value(node_index_000, n_components);
375 auto value_100 = get_value(node_index_100, n_components);
376
377 auto value_010 = get_value(node_index_010, n_components);
378 auto value_110 = get_value(node_index_110, n_components);
379
380 auto value_001 = get_value(node_index_001, n_components);
381 auto value_101 = get_value(node_index_101, n_components);
382
383 auto value_011 = get_value(node_index_011, n_components);
384 auto value_111 = get_value(node_index_111, n_components);
385
386 // Interpolate
387 for (unsigned int c : std::views::iota(0U, n_components))
388 {
389 // cppcheck-suppress-begin containerOutOfBounds
390 value[c] =
391 (1.0 - weights[0]) * (1.0 - weights[1]) * (1.0 - weights[2]) * value_000[c] +
392 weights[0] * (1.0 - weights[1]) * (1.0 - weights[2]) * value_100[c] +
393 (1.0 - weights[0]) * weights[1] * (1.0 - weights[2]) * value_010[c] +
394 weights[0] * weights[1] * (1.0 - weights[2]) * value_110[c] +
395 (1.0 - weights[0]) * (1.0 - weights[1]) * weights[2] * value_001[c] +
396 weights[0] * (1.0 - weights[1]) * weights[2] * value_101[c] +
397 (1.0 - weights[0]) * weights[1] * weights[2] * value_011[c] +
398 weights[0] * weights[1] * weights[2] * value_111[c];
399 // cppcheck-suppress-end containerOutOfBounds
400 }
401 }
402
403 return value;
404}
405
406template <unsigned int dim, typename number>
407inline number
408ReadBinary<dim, number>::get_scalar_value(const dealii::Point<dim> &point,
409 [[maybe_unused]] const std::string &scalar_name)
410{
411 Assert(n_values == n_points,
412 dealii::ExcMessage("The number of points should match the number of values in a "
413 "binary file for a scalar field. Make sure the file size is "
414 "correct and you are trying to access a scalar field."));
415 return interpolate(point, 1)[0];
416}
417
418template <unsigned int dim, typename number>
419inline dealii::Vector<number>
420ReadBinary<dim, number>::get_vector_value(const dealii::Point<dim> &point,
421 [[maybe_unused]] const std::string &vector_name)
422{
423 Assert((n_values / dim) == n_points,
424 dealii::ExcMessage("The number of points should match the number of values "
425 "divided by the dimension in a binary file for a vector "
426 "field. Make sure the file size is correct and you are "
427 "trying to access a vector field."));
428
429 return interpolate(point, dim);
430}
431
432template <unsigned int dim, typename number>
433inline void
435{
436 for (dealii::types::global_dof_index i : std::views::iota(0U, n_values))
437 {
438 Assert(i < data.size(),
439 dealii::ExcMessage("Index out of bounds in ReadBinary::print_file"));
440 ConditionalOStreams::pout_summary() << this->data.at(i) << "\n";
441 }
442 ConditionalOStreams::pout_summary() << std::flush;
443}
444
445PRISMS_PF_END_NAMESPACE
static dealii::ConditionalOStream & pout_summary()
Log output stream for writing a summary.log file.
Definition conditional_ostreams.cc:35
dealii::Vector< number > get_value(const dealii::types::global_dof_index index, const unsigned int n_components)
Get vector value for a given index.
Definition read_binary.h:204
std::vector< number > data
Data array to hold the read in values.
Definition read_binary.h:109
ReadBinary(const InitialConditionFile &_ic_file, const SpatialDiscretization< dim > &_spatial_discretization)
Constructor.
Definition read_binary.h:113
number get_scalar_value(const dealii::Point< dim > &point, const std::string &scalar_name) override
Get scalar value for a given point.
Definition read_binary.h:408
void check_file_size()
Check the size of the binary file and make sure it matches the expected size (in bytes).
Definition read_binary.h:157
dealii::Vector< number > interpolate(const dealii::Point< dim > &point, const unsigned int n_components)
Get vector value for a given point.
Definition read_binary.h:228
dealii::types::global_dof_index n_values
Number of values (n_points * n_components).
Definition read_binary.h:104
dealii::types::global_dof_index n_points
Number of grid points.
Definition read_binary.h:99
static void write_file(const std::vector< number > &data, const InitialConditionFile &ic_file)
Write a binary file for testing.
Definition read_binary.h:188
void print_file() override
Print the binary file to text for debugging.
Definition read_binary.h:434
dealii::Vector< number > get_vector_value(const dealii::Point< dim > &point, const std::string &vector_name) override
Get vector value for a given point.
Definition read_binary.h:420
const InitialConditionFile & ic_file
Initial condition file object.
Definition read_field_base.h:86
ReadFieldBase(const InitialConditionFile &_ic_file, const SpatialDiscretization< dim > &_spatial_discretization)
Constructor.
Definition read_field_base.h:90
const SpatialDiscretization< dim > & spatial_discretization
Spatial discretization object.
Definition read_field_base.h:81
Definition conditional_ostreams.cc:20
@ ReadBinary
Definition read_field_factory.h:32
@ Rectangular
Definition spatial_discretization.h:27
Struct that store the read-in information for a single file.
Definition load_initial_condition_parameters.h:25
Struct that holds spatial discretization parameters.
Definition spatial_discretization.h:230
TriangulationType type
Definition spatial_discretization.h:291
@ FlatBinary
Definition type_enums.h:134