Loading...
Searching...
No Matches
radial_profiles.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <numeric>
4
8#include "util/table.hpp"
9
10namespace terra::shell {
11
12/// @brief Simple struct, holding device-views for radial profiles.
13///
14/// See `radial_profiles()` and `radial_profiles_to_table()` for radial profiles computation.
15template < typename ScalarType >
17{
18 explicit RadialProfiles( int radial_shells )
19 : radial_min_( "radial_profiles_min", radial_shells )
20 , radial_max_( "radial_profiles_max", radial_shells )
21 , radial_sum_( "radial_profiles_sum", radial_shells )
22 , radial_cnt_( "radial_profiles_cnt", radial_shells )
23 , radial_avg_( "radial_profiles_avg", radial_shells )
24 {}
25
31};
32
33/// @brief Compute radial profiles (min, max, sum, count, average) for a Q1 scalar field (on device!).
34///
35/// @details Computes the radial profiles for a Q1 scalar field by iterating over all radial shells.
36/// The profiles include minimum, maximum, sum, and count of nodes in each shell.
37/// The radial shells are defined by the radial dimension of the Q1 scalar field.
38/// The output is a RadialProfiles struct with device-side arrays of size num_shells which contain:
39/// - Minimum value in the shell
40/// - Maximum value in the shell
41/// - Sum of values in the shell
42/// - Count of nodes in the shell
43/// - Average of nodes in the shell (sum / count)
44/// Performs reduction per process on the device and also inter-device reduction with MPI_Allreduce.
45/// All processes will carry the same result.
46///
47/// Mask data is used internally to filter out non-owned nodes. So this really only reduces owned nodes.
48///
49/// @note The returned Views are still on the device.
50/// To convert it to a util::Table for output, use the `radial_profiles_to_table()` function.
51/// So unless you really want to compute any further data on the device using the profile, you like always want
52/// to call the `radial_profiles_to_table()` function as outlined below.
53///
54/// To nicely format the output, use the `radial_profiles_to_table()` function, e.g., via
55/// @code
56/// auto radii = domain_info.radii(); // the DomainInfo is also available in the DistributedDomain
57/// auto table = radial_profiles_to_table( radial_profiles( some_field_like_temperature ), radii );
58/// std::ofstream out( "radial_profiles.csv" );
59/// table.print_csv( out );
60/// @endcode
61///
62/// @tparam ScalarType Scalar type of the field.
63/// @param data Q1 scalar field data.
64/// @param subdomain_shell_idx global shell indices array with layout
65/// @code subdomain_shell_idx( local_subdomain_id, node_r_idx ) = global_shell_idx @endcode
66/// compute with \ref terra::grid::shell::subdomain_shell_idx.
67/// @param num_global_shells number of global shells
68/// @return RadialProfiles struct containing [min, max, sum, count, average] for each radial shell (still on the device).
69template < typename ScalarType >
72 const grid::Grid2DDataScalar< int >& subdomain_shell_idx,
73 const int num_global_shells )
74{
76
77 const auto data_grid = data.grid_data();
78 const auto data_mask = data.mask_data();
79
80 if ( data_grid.extent( 0 ) != subdomain_shell_idx.extent( 0 ) ||
81 data_grid.extent( 3 ) != subdomain_shell_idx.extent( 1 ) )
82 {
83 Kokkos::abort( "radial_profiles: Data and subdomain_shell_idx do not have matching dimensions." );
84 }
85
86 Kokkos::parallel_for(
87 "radial profiles init", num_global_shells, KOKKOS_LAMBDA( int r ) {
88 radial_profiles.radial_min_( r ) = Kokkos::Experimental::finite_max_v< ScalarType >;
89 radial_profiles.radial_max_( r ) = Kokkos::Experimental::finite_min_v< ScalarType >;
90 radial_profiles.radial_sum_( r ) = 0;
91 radial_profiles.radial_cnt_( r ) = 0;
92 } );
93
94 Kokkos::fence();
95
96 Kokkos::parallel_for(
97 "radial profiles reduction",
98 Kokkos::MDRangePolicy(
99 { 0, 0, 0, 0 },
100 { data_grid.extent( 0 ), data_grid.extent( 1 ), data_grid.extent( 2 ), data_grid.extent( 3 ) } ),
101 KOKKOS_LAMBDA( int local_subdomain_id, int x, int y, int r ) {
102 if ( !util::has_flag( data_mask( local_subdomain_id, x, y, r ), grid::NodeOwnershipFlag::OWNED ) )
103 {
104 return;
105 }
106
107 const int global_shell_idx = subdomain_shell_idx( local_subdomain_id, r );
108
109 Kokkos::atomic_min(
110 &radial_profiles.radial_min_( global_shell_idx ), data_grid( local_subdomain_id, x, y, r ) );
111 Kokkos::atomic_max(
112 &radial_profiles.radial_max_( global_shell_idx ), data_grid( local_subdomain_id, x, y, r ) );
113 Kokkos::atomic_add(
114 &radial_profiles.radial_sum_( global_shell_idx ), data_grid( local_subdomain_id, x, y, r ) );
115 Kokkos::atomic_add( &radial_profiles.radial_cnt_( global_shell_idx ), static_cast< ScalarType >( 1 ) );
116 } );
117
118 Kokkos::fence();
119
120 // Perform MPI_Allreduce on host to avoid GPU-aware MPI issues with MPI_IN_PLACE
121 // at large rank counts. The buffers are small (num_global_shells doubles), so the
122 // device-to-host copy overhead is negligible.
123 {
124 auto host_min = Kokkos::create_mirror_view_and_copy( Kokkos::HostSpace(), radial_profiles.radial_min_ );
125 auto host_max = Kokkos::create_mirror_view_and_copy( Kokkos::HostSpace(), radial_profiles.radial_max_ );
126 auto host_sum = Kokkos::create_mirror_view_and_copy( Kokkos::HostSpace(), radial_profiles.radial_sum_ );
127 auto host_cnt = Kokkos::create_mirror_view_and_copy( Kokkos::HostSpace(), radial_profiles.radial_cnt_ );
128
129 MPI_Allreduce(
130 MPI_IN_PLACE, host_min.data(), host_min.size(), mpi::mpi_datatype< ScalarType >(), MPI_MIN, MPI_COMM_WORLD );
131 MPI_Allreduce(
132 MPI_IN_PLACE, host_max.data(), host_max.size(), mpi::mpi_datatype< ScalarType >(), MPI_MAX, MPI_COMM_WORLD );
133 MPI_Allreduce(
134 MPI_IN_PLACE, host_sum.data(), host_sum.size(), mpi::mpi_datatype< ScalarType >(), MPI_SUM, MPI_COMM_WORLD );
135 MPI_Allreduce(
136 MPI_IN_PLACE, host_cnt.data(), host_cnt.size(), mpi::mpi_datatype< ScalarType >(), MPI_SUM, MPI_COMM_WORLD );
137
138 Kokkos::deep_copy( radial_profiles.radial_min_, host_min );
139 Kokkos::deep_copy( radial_profiles.radial_max_, host_max );
140 Kokkos::deep_copy( radial_profiles.radial_sum_, host_sum );
141 Kokkos::deep_copy( radial_profiles.radial_cnt_, host_cnt );
142 }
143
144 Kokkos::parallel_for(
145 "radial profiles avg", num_global_shells, KOKKOS_LAMBDA( int r ) {
146 radial_profiles.radial_avg_( r ) = radial_profiles.radial_sum_( r ) / radial_profiles.radial_cnt_( r );
147 } );
148
149 Kokkos::fence();
150
151 return radial_profiles;
152}
153
154/// @brief Convert radial profile data to a util::Table for analysis or output.
155///
156/// @details Converts the radial profile data (min, max, sum, avg, count) into a util::Table.
157///
158/// This table can then be used for further analysis or output to CSV/JSON.
159/// The table will have the following columns:
160/// - tag: "radial_profiles"
161/// - shell_idx: Index of the radial shell
162/// - radius: Radius of the shell
163/// - min: Minimum value in the shell
164/// - max: Maximum value in the shell
165/// - sum: Average value in the shell
166/// - avg: Average value in the shell
167/// - cnt: Count of nodes in the shell
168///
169/// To use this function, you can compute the radial profiles using `radial_profiles()` and then convert
170/// the result to a table using this function:
171///
172/// @code
173/// auto radii = domain_info.radii(); // the DomainInfo is also available in the DistributedDomain
174/// auto table = radial_profiles_to_table( radial_profiles( some_field_like_temperature ), radii );
175/// std::ofstream out( "radial_profiles.csv" );
176/// table.print_csv( out );
177/// @endcode
178///
179/// @tparam ScalarType Scalar type of the profile data.
180/// @param radial_profiles RadialProfiles struct containing radial profile statistics.
181/// Compute this using `radial_profiles()` function. Data is expected to be on the device still.
182/// It is copied to host for table creation in this function.
183/// @param radii Vector of shell radii. Can for instance be obtained from the DomainInfo.
184/// @return Table with columns: tag, shell_idx, radius, min, max, avg, cnt.
185template < typename ScalarType >
187 radial_profiles_to_table( const RadialProfiles< ScalarType >& radial_profiles, const std::vector< double >& radii )
188{
189 if ( radii.size() != radial_profiles.radial_min_.extent( 0 ) ||
190 radii.size() != radial_profiles.radial_max_.extent( 0 ) ||
191 radii.size() != radial_profiles.radial_sum_.extent( 0 ) ||
192 radii.size() != radial_profiles.radial_cnt_.extent( 0 ) )
193 {
194 throw std::runtime_error( "Radial profiles and radii do not have the same number of shells." );
195 }
196
197 const auto radial_profiles_host_min =
198 Kokkos::create_mirror_view_and_copy( Kokkos::HostSpace(), radial_profiles.radial_min_ );
199 const auto radial_profiles_host_max =
200 Kokkos::create_mirror_view_and_copy( Kokkos::HostSpace(), radial_profiles.radial_max_ );
201 const auto radial_profiles_host_sum =
202 Kokkos::create_mirror_view_and_copy( Kokkos::HostSpace(), radial_profiles.radial_sum_ );
203 const auto radial_profiles_host_cnt =
204 Kokkos::create_mirror_view_and_copy( Kokkos::HostSpace(), radial_profiles.radial_cnt_ );
205 const auto radial_profiles_host_avg =
206 Kokkos::create_mirror_view_and_copy( Kokkos::HostSpace(), radial_profiles.radial_avg_ );
207
208 util::Table table;
209 for ( int r = 0; r < radii.size(); r++ )
210 {
211 table.add_row(
212 { { "tag", "radial_profiles" },
213 { "shell_idx", r },
214 { "radius", radii[r] },
215 { "min", radial_profiles_host_min( r ) },
216 { "max", radial_profiles_host_max( r ) },
217 { "sum", radial_profiles_host_sum( r ) },
218 { "avg", radial_profiles_host_avg( r ) },
219 { "cnt", radial_profiles_host_cnt( r ) } } );
220 }
221 return table;
222}
223
224/// @brief Linearly interpolates radial data from a std::vector (host) into a 2D grid (device) with dimensions
225/// (local_subdomain_id, subdomain_size_r) for straightforward use in kernels.
226///
227/// This function takes care of handling possibly duplicated nodes at subdomain boundaries.
228///
229/// @note This will clamp values outside the passed radial profile values to the first and last values in the vector
230/// respectively.
231///
232/// @param profile_out_label label of the returned Kokkos::View
233/// @param coords_radii radii of each node for all local subdomains - see \ref terra::grid::shell::subdomain_shell_radii(),
234/// the output grid data will have the same extents
235/// @param profile_radii input profile: vector of radii
236/// @param profile_values input profile: vector of values to interpolate
237/// @return Kokkos::View with the same dimensions of coords_radii, populated with interpolated values per subdomain
238template <
239 std::floating_point GridDataType,
240 std::floating_point ProfileInRadiiType,
241 std::floating_point ProfileInValueType,
242 std::floating_point ProfileOutDataType = double >
244 const std::string& profile_out_label,
245 const grid::Grid2DDataScalar< GridDataType >& coords_radii,
246 const std::vector< ProfileInRadiiType >& profile_radii,
247 const std::vector< ProfileInValueType >& profile_values )
248{
250 profile_out_label, coords_radii.extent( 0 ), coords_radii.extent( 1 ) );
251 auto profile_data_host = Kokkos::create_mirror_view( Kokkos::HostSpace{}, profile_data );
252
253 auto coords_radii_host = Kokkos::create_mirror_view_and_copy( Kokkos::HostSpace{}, coords_radii );
254
255 // Sort ascending (by radius - from innermost to outermost shell).
256
257 std::vector< int > idx( profile_radii.size() );
258 std::iota( idx.begin(), idx.end(), 0 );
259
260 // sort indices by comparing keys
261 std::sort( idx.begin(), idx.end(), [&]( int a, int b ) { return profile_radii[a] < profile_radii[b]; } );
262
263 // apply permutation
264 std::vector< ProfileInRadiiType > profile_radii_sorted( profile_radii );
265 std::vector< ProfileInValueType > profile_values_sorted( profile_values );
266
267 for ( size_t i = 0; i < idx.size(); i++ )
268 {
269 profile_radii_sorted[i] = profile_radii[idx[i]];
270 profile_values_sorted[i] = profile_values[idx[i]];
271 }
272
273 for ( int local_subdomain_id = 0; local_subdomain_id < profile_data_host.extent( 0 ); local_subdomain_id++ )
274 {
275 const size_t subdomain_size_r = profile_data_host.extent( 1 );
276 const size_t input_profile_size = profile_radii_sorted.size();
277
278 if ( input_profile_size == 0 )
279 {
280 continue;
281 }
282
283 if ( input_profile_size == 1 )
284 {
285 for ( int r = 0; r < subdomain_size_r; r++ )
286 {
287 profile_data_host( local_subdomain_id, r ) = profile_values_sorted[0];
288 }
289 continue;
290 }
291
292 size_t i = 0;
293
294 for ( size_t r = 0; r < subdomain_size_r; ++r )
295 {
296 double x = coords_radii_host( local_subdomain_id, r );
297
298 // Advance src index to the correct interval
299 while ( i + 1 < input_profile_size && profile_radii_sorted[i + 1] < x )
300 {
301 ++i;
302 }
303
304 // If x is below the first point → clamp
305 if ( x <= profile_radii_sorted.front() )
306 {
307 profile_data_host( local_subdomain_id, r ) = profile_values_sorted.front();
308 continue;
309 }
310
311 // If x is above the last point → clamp
312 if ( x >= profile_radii_sorted.back() )
313 {
314 profile_data_host( local_subdomain_id, r ) = profile_values_sorted.back();
315 continue;
316 }
317
318 // Interpolate between i and i+1
319 profile_data_host( local_subdomain_id, r ) = util::interpolate_linear_1D(
320 profile_radii_sorted[i],
321 profile_radii_sorted[i + 1],
322 profile_values_sorted[i],
323 profile_values_sorted[i + 1],
324 x,
325 true );
326 }
327 }
328
329 Kokkos::deep_copy( profile_data, profile_data_host );
330
331 return profile_data;
332}
333
334/// @brief Interpolating a radial profile from a CSV file into a grid.
335///
336/// This is just a convenient wrapper around \ref interpolate_radial_profile_into_subdomains() that also
337/// reads the CSV file.
338///
339/// @note This will clamp values outside the passed radial profile values to the first and last values in the vector
340/// respectively.
341///
342/// @param filename csv file name containing the radial profile (radii and at least one value column)
343/// @param key_radii name of the radii column
344/// @param key_values name of the value column
345/// @param coords_radii radii of each node for all local subdomains - see \ref terra::grid::shell::subdomain_shell_radii(),
346/// the output grid data will have the same extents
347/// @return Kokkos::View with the same dimensions of coords_radii, populated with interpolated values per subdomain
348template < std::floating_point GridDataType, std::floating_point ProfileOutDataType = double >
350 const std::string& filename,
351 const std::string& key_radii,
352 const std::string& key_values,
353 const grid::Grid2DDataScalar< GridDataType >& coords_radii )
354{
355 auto profile_table_result = util::read_table_from_csv( filename );
356 if ( profile_table_result.is_err() )
357 {
358 util::logroot << profile_table_result.error() << std::endl;
359 Kokkos::abort( "Error reading csv file into table." );
360 }
361 const auto& profile_table = profile_table_result.unwrap();
362
363 const auto profile_radii = profile_table.column_as_vector< double >( key_radii );
364 const auto profile_values = profile_table.column_as_vector< double >( key_values );
365
367 "radial_profile_" + key_values, coords_radii, profile_radii, profile_values );
368}
369
370/// @brief Interpolating a constant radial profile into a grid.
371///
372/// @param coords_radii radii of each node for all local subdomains - see \ref terra::grid::shell::subdomain_shell_radii(),
373/// the output grid data will have the same extents
374/// @param value the constant value to interpolate
375/// @return Kokkos::View with the same dimensions of coords_radii, populated with interpolated values per subdomain
376template < std::floating_point GridDataType, std::floating_point ProfileOutDataType = double >
378 const grid::Grid2DDataScalar< GridDataType >& coords_radii,
379 const ProfileOutDataType& value )
380{
382 "radial_profile_constant", coords_radii.extent( 0 ), coords_radii.extent( 1 ) );
383 kernels::common::set_constant( profile_data, value );
384 return profile_data;
385}
386
387} // namespace terra::shell
const grid::Grid4DDataScalar< ScalarType > & grid_data() const
Get const reference to grid data.
Definition vector_q1.hpp:139
const grid::Grid4DDataScalar< grid::NodeOwnershipFlag > & mask_data() const
Get const reference to mask data.
Definition vector_q1.hpp:144
Table class for storing and manipulating tabular data.
Definition table.hpp:84
void add_row(const Row &row_data)
Add a row to the table. Adds "id" and "timestamp" columns automatically.
Definition table.hpp:140
Kokkos::View< ScalarType *, Layout > Grid1DDataScalar
Definition grid_types.hpp:18
Kokkos::View< ScalarType **, Layout > Grid2DDataScalar
Definition grid_types.hpp:21
void set_constant(const grid::Grid2DDataScalar< ScalarType > &x, ScalarType value)
Definition grid_operations.hpp:12
Definition radial_profiles.hpp:10
grid::Grid2DDataScalar< ProfileOutDataType > interpolate_constant_radial_profile(const grid::Grid2DDataScalar< GridDataType > &coords_radii, const ProfileOutDataType &value)
Interpolating a constant radial profile into a grid.
Definition radial_profiles.hpp:377
util::Table radial_profiles_to_table(const RadialProfiles< ScalarType > &radial_profiles, const std::vector< double > &radii)
Convert radial profile data to a util::Table for analysis or output.
Definition radial_profiles.hpp:187
grid::Grid2DDataScalar< ProfileOutDataType > interpolate_radial_profile_into_subdomains_from_csv(const std::string &filename, const std::string &key_radii, const std::string &key_values, const grid::Grid2DDataScalar< GridDataType > &coords_radii)
Interpolating a radial profile from a CSV file into a grid.
Definition radial_profiles.hpp:349
RadialProfiles< ScalarType > radial_profiles(const linalg::VectorQ1Scalar< ScalarType > &data, const grid::Grid2DDataScalar< int > &subdomain_shell_idx, const int num_global_shells)
Compute radial profiles (min, max, sum, count, average) for a Q1 scalar field (on device!...
Definition radial_profiles.hpp:70
grid::Grid2DDataScalar< ProfileOutDataType > interpolate_radial_profile_into_subdomains(const std::string &profile_out_label, const grid::Grid2DDataScalar< GridDataType > &coords_radii, const std::vector< ProfileInRadiiType > &profile_radii, const std::vector< ProfileInValueType > &profile_values)
Linearly interpolates radial data from a std::vector (host) into a 2D grid (device) with dimensions (...
Definition radial_profiles.hpp:243
Result< Table > read_table_from_csv(const std::string &filename)
Attempts to read a csv file and converts that into a Table instance.
Definition table.hpp:523
double interpolate_linear_1D(const double pos_a, const double pos_b, const double val_a, const double val_b, const double pos, const bool clamp)
Computes the linearly interpolated value at a specified point using two surrounding data points.
Definition interpolation.hpp:17
constexpr bool has_flag(E mask_value, E flag) noexcept
Checks if a bitmask value contains a specific flag.
Definition bit_masking.hpp:43
detail::PrefixCout logroot([]() { return detail::log_prefix();})
std::ostream subclass that just logs on root and adds a timestamp for each line.
Simple struct, holding device-views for radial profiles.
Definition radial_profiles.hpp:17
grid::Grid1DDataScalar< ScalarType > radial_cnt_
Definition radial_profiles.hpp:29
grid::Grid1DDataScalar< ScalarType > radial_sum_
Definition radial_profiles.hpp:28
grid::Grid1DDataScalar< ScalarType > radial_avg_
Definition radial_profiles.hpp:30
grid::Grid1DDataScalar< ScalarType > radial_min_
Definition radial_profiles.hpp:26
RadialProfiles(int radial_shells)
Definition radial_profiles.hpp:18
grid::Grid1DDataScalar< ScalarType > radial_max_
Definition radial_profiles.hpp:27