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:
39
43 void
44 print_file() override;
45
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>
117{
118 // Make sure the dataset format is correct
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
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 {
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);
144 dealii::ExcMessage("Could not open binary file: " +
145 this->ic_file.filename));
146
147 // Reserve space in the data vector
148 data.reserve(n_values);
149
150 // Read in the data
151 for (dealii::types::global_dof_index i : std::views::iota(0U, n_values))
152 {
153 std::array<char, sizeof(number)> buffer;
154 data_file.read(buffer.data(), sizeof(number));
155 data.push_back(std::bit_cast<number>(buffer));
156 }
157 data_file.close();
158}
159
160template <unsigned int dim, typename number>
161inline void
163{
164 // Grab the file size of the binary file in bytes
165 auto file_size = std::filesystem::file_size(this->ic_file.filename);
166
167 // Compute the expected size of the binary file. This is simply the number of points
168 // multiplied by the size of each point in bytes.
169 auto expected_size_scalar = static_cast<std::uintmax_t>(n_points * sizeof(number));
170 auto expected_size_vector = static_cast<std::uintmax_t>(dim * expected_size_scalar);
171
172 // Make sure expected size is not zero
175 dealii::ExcMessage(
176 "Expected input array size is zero, check that the number of data points "
177 "in each used direction is set correctly in the input file for your binary file. "
178 "You likely have the number of data points set to zero in all directions."));
179 // Make sure the size matches for either a scalar or vector
181 dealii::ExcMessage(
182 "Expected binary file size (" + std::to_string(expected_size_scalar) +
183 " bytes for scalar or " + std::to_string(expected_size_vector) +
184 " bytes for vector) does not match actual file size (" +
185 std::to_string(file_size) + " bytes)."));
186
187 // Set the number of values
188 n_values = file_size == expected_size_scalar ? n_points : dim * n_points;
189}
190
191template <unsigned int dim, typename number>
192inline void
193ReadBinary<dim, number>::write_file(const std::vector<number> &data,
194 const InitialConditionFile &ic_file)
195{
196 // Try to open the file
197 std::ofstream data_file(ic_file.filename, std::ios::binary);
199 dealii::ExcMessage("Could not open binary file: " + ic_file.filename));
200
201 // Write the data
202 for (dealii::types::global_dof_index j : std::views::iota(0U, data.size()))
203 {
204 auto buffer = std::bit_cast<std::array<char, sizeof(number)>>(data[j]);
205 data_file.write(buffer.data(), sizeof(number));
206 }
207 data_file.close();
208}
209
210template <unsigned int dim, typename number>
211inline dealii::Vector<number>
212ReadBinary<dim, number>::get_value(const dealii::types::global_dof_index index,
213 const unsigned int n_components)
214{
215 // Create vector to hold the value with the correct number of components
216 dealii::Vector<number> value(n_components);
217
218 // Check that the number of components matches what we expect from the number of values
219 Assert(n_values % n_components == 0,
220 dealii::ExcMessage("The number of components requested in the get_value call "
221 "does not match the number of values in the binary file."));
222
223 // Fill the value vector
224 for (unsigned int component : std::views::iota(0U, n_components))
225 {
226 Assert(((n_components * index) + component) < data.size(),
227 dealii::ExcMessage("Index out of bounds in ReadBinary::get_vals"));
228 value[component] = data[(n_components * index) + component];
229 }
230
231 return value;
232}
233
234template <unsigned int dim, typename number>
235inline dealii::Vector<number>
237 const unsigned int n_components)
238{
240 dealii::ExcMessage(
241 "Number of components for interpolation must be 1 (scalar) or dim (vector)"));
242
243 // Create a vector to hold the interpolated value
244 dealii::Vector<number> value(n_components);
245
246 // Get the spatial discretization
247 const auto &spatial_discretization = this->spatial_discretization;
248
249 // Compute the spacing in each direction
250 std::array<number, dim> spacing;
251 for (unsigned int d : std::views::iota(0U, dim))
252 {
253 spacing[d] = spatial_discretization.rectangular_mesh.size[d] /
254 static_cast<number>(this->ic_file.n_data_points[d]);
255 }
256
257 // Compute the indices of the lower corner of the cell containing the point
258 std::array<dealii::types::global_dof_index, dim> lower_indices;
259 for (unsigned int d : std::views::iota(0U, dim))
260 {
262 static_cast<dealii::types::global_dof_index>(std::floor(point[d] / spacing[d]));
263 // Make sure we don't go out of bounds
264 if (lower_indices[d] >= this->ic_file.n_data_points[d] - 1)
265 {
266 lower_indices[d] = this->ic_file.n_data_points[d] - 2;
267 }
268 }
269
270 // Compute the weights for interpolation in each direction
271 std::array<number, dim> weights;
272 for (unsigned int d : std::views::iota(0U, dim))
273 {
274 weights[d] = (point[d] - lower_indices[d] * spacing[d]) / spacing[d];
275 }
276
277 // Perform multilinear interpolation based on the dimension
278 if constexpr (dim == 1)
279 {
280 // Here is the map of the nodes in 1D:
281 // 0 — 1
282
283 // Grab the index of the lower left corner of the cell
284 // cppcheck-suppress-begin containerOutOfBounds
285 auto lower_index = lower_indices[0];
286 // cppcheck-suppress-end containerOutOfBounds
287
288 // Get the indices of the two points in 1D
290 auto node_index_1 = lower_index + 1;
291
292 // Get the values at these points
293 auto value_0 = get_value(node_index_0, n_components);
294 auto value_1 = get_value(node_index_1, n_components);
295
296 // Interpolate
297 for (unsigned int c : std::views::iota(0U, n_components))
298 {
299 // cppcheck-suppress-begin containerOutOfBounds
300 value[c] = (1.0 - weights[0]) * value_0[c] + weights[0] * value_1[c];
301 // cppcheck-suppress-end containerOutOfBounds
302 }
303 }
304 else if constexpr (dim == 2)
305 {
306 // Here is the map of the nodes in 2D:
307 // 01 — 11
308 // | |
309 // 00 — 10
310
311 // Grab the row length in the x-direction (0th direction)
312 auto row_length_0 = this->ic_file.n_data_points[0];
313
314 // Grab the index of the lower left corner of the cell
315 // cppcheck-suppress-begin containerOutOfBounds
317 // cppcheck-suppress-end containerOutOfBounds
318
319 // Get the indices of the four points in 2D
321 auto node_index_10 = lower_index + 1;
322
325
326 // Get the values at these points
327 auto value_00 = get_value(node_index_00, n_components);
328 auto value_10 = get_value(node_index_10, n_components);
329
330 auto value_01 = get_value(node_index_01, n_components);
331 auto value_11 = get_value(node_index_11, n_components);
332
333 // Interpolate
334 for (unsigned int c : std::views::iota(0U, n_components))
335 {
336 // cppcheck-suppress-begin containerOutOfBounds
337 value[c] = (1.0 - weights[0]) * (1.0 - weights[1]) * value_00[c] +
338 weights[0] * (1.0 - weights[1]) * value_10[c] +
339 (1.0 - weights[0]) * weights[1] * value_01[c] +
340 weights[0] * weights[1] * value_11[c];
341 // cppcheck-suppress-end containerOutOfBounds
342 }
343 }
344 else if constexpr (dim == 3)
345 {
346 // Here is the map of the nodes in 3D:
347 //
348 // 011 ———— 111
349 // / | / |
350 // / 001 / |
351 // 010 ——— 110 101
352 // | / | /
353 // | / | /
354 // 000 ———— 100
355
356 // Grab the row length in the x-direction (0th direction)
357 auto row_length_0 = this->ic_file.n_data_points[0];
358 // Grab the row length in the y-direction (1st direction)
359 auto row_length_1 = this->ic_file.n_data_points[1];
360
361 // Grab the index of the lower left corner of the cell
362 // cppcheck-suppress-begin containerOutOfBounds
365 // cppcheck-suppress-end containerOutOfBounds
366
367 // Get the indices of the eight points in 3D
369 auto node_index_100 = lower_index + 1;
370
373
376
378 auto node_index_111 =
380
381 // Get the values at these points
382 auto value_000 = get_value(node_index_000, n_components);
383 auto value_100 = get_value(node_index_100, n_components);
384
385 auto value_010 = get_value(node_index_010, n_components);
386 auto value_110 = get_value(node_index_110, n_components);
387
388 auto value_001 = get_value(node_index_001, n_components);
389 auto value_101 = get_value(node_index_101, n_components);
390
391 auto value_011 = get_value(node_index_011, n_components);
392 auto value_111 = get_value(node_index_111, n_components);
393
394 // Interpolate
395 for (unsigned int c : std::views::iota(0U, n_components))
396 {
397 // cppcheck-suppress-begin containerOutOfBounds
398 value[c] =
399 (1.0 - weights[0]) * (1.0 - weights[1]) * (1.0 - weights[2]) * value_000[c] +
400 weights[0] * (1.0 - weights[1]) * (1.0 - weights[2]) * value_100[c] +
401 (1.0 - weights[0]) * weights[1] * (1.0 - weights[2]) * value_010[c] +
402 weights[0] * weights[1] * (1.0 - weights[2]) * value_110[c] +
403 (1.0 - weights[0]) * (1.0 - weights[1]) * weights[2] * value_001[c] +
404 weights[0] * (1.0 - weights[1]) * weights[2] * value_101[c] +
405 (1.0 - weights[0]) * weights[1] * weights[2] * value_011[c] +
406 weights[0] * weights[1] * weights[2] * value_111[c];
407 // cppcheck-suppress-end containerOutOfBounds
408 }
409 }
410
411 return value;
412}
413
414template <unsigned int dim, typename number>
415inline number
417 [[maybe_unused]] const std::string &scalar_name)
418{
419 Assert(n_values == n_points,
420 dealii::ExcMessage("The number of points should match the number of values in a "
421 "binary file for a scalar field. Make sure the file size is "
422 "correct and you are trying to access a scalar field."));
423 return interpolate(point, 1)[0];
424}
425
426template <unsigned int dim, typename number>
427inline dealii::Vector<number>
429 [[maybe_unused]] const std::string &vector_name)
430{
431 Assert((n_values / dim) == n_points,
432 dealii::ExcMessage("The number of points should match the number of values "
433 "divided by the dimension in a binary file for a vector "
434 "field. Make sure the file size is correct and you are "
435 "trying to access a vector field."));
436
437 return interpolate(point, dim);
438}
439
440template <unsigned int dim, typename number>
441inline void
443{
444 for (dealii::types::global_dof_index i : std::views::iota(0U, n_values))
445 {
446 Assert(i < data.size(),
447 dealii::ExcMessage("Index out of bounds in ReadBinary::print_file"));
448 ConditionalOStreams::pout_summary() << this->data.at(i) << "\n";
449 }
450 ConditionalOStreams::pout_summary() << std::flush;
451}
452
static dealii::ConditionalOStream & pout_summary()
Log output stream for writing a summary.log file.
Definition conditional_ostreams.cc:34
Class to read in a flat binary file and provide values at given points.
Definition read_binary.h:32
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:212
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:416
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:162
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:236
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:193
void print_file() override
Print the binary file to text for debugging.
Definition read_binary.h:442
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:428
Definition read_field_base.h:23
const InitialConditionFile & ic_file
Initial condition file object.
Definition read_field_base.h:86
@ Value
Use value of the variable as a criterion for refinement.
Definition grid_refiner_criterion.h:31
Definition conditional_ostreams.cc:20
@ Rectangular
Definition spatial_discretization.h:27
Struct that store the read-in information for a single file.
Definition load_initial_condition_parameters.h:25
std::array< dealii::types::global_dof_index, 3 > n_data_points
Definition load_initial_condition_parameters.h:39
std::vector< std::string > file_variable_names
Definition load_initial_condition_parameters.h:33
DataFormatType dataset_format
Definition load_initial_condition_parameters.h:30
std::string filename
Definition load_initial_condition_parameters.h:27
Struct that holds spatial discretization parameters.
Definition spatial_discretization.h:230
@ FlatBinary
Definition type_enums.h:112