15template <
typename ScalarType >
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 )
69template <
typename ScalarType >
73 const int num_global_shells )
80 if ( data_grid.extent( 0 ) != subdomain_shell_idx.extent( 0 ) ||
81 data_grid.extent( 3 ) != subdomain_shell_idx.extent( 1 ) )
83 Kokkos::abort(
"radial_profiles: Data and subdomain_shell_idx do not have matching dimensions." );
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 >;
97 "radial profiles reduction",
98 Kokkos::MDRangePolicy(
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 ) {
107 const int global_shell_idx = subdomain_shell_idx( local_subdomain_id, r );
110 &
radial_profiles.radial_min_( global_shell_idx ), data_grid( local_subdomain_id, x, y, r ) );
112 &
radial_profiles.radial_max_( global_shell_idx ), data_grid( local_subdomain_id, x, y, r ) );
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 ) );
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_ );
130 MPI_IN_PLACE, host_min.data(), host_min.size(), mpi::mpi_datatype< ScalarType >(), MPI_MIN, MPI_COMM_WORLD );
132 MPI_IN_PLACE, host_max.data(), host_max.size(), mpi::mpi_datatype< ScalarType >(), MPI_MAX, MPI_COMM_WORLD );
134 MPI_IN_PLACE, host_sum.data(), host_sum.size(), mpi::mpi_datatype< ScalarType >(), MPI_SUM, MPI_COMM_WORLD );
136 MPI_IN_PLACE, host_cnt.data(), host_cnt.size(), mpi::mpi_datatype< ScalarType >(), MPI_SUM, MPI_COMM_WORLD );
144 Kokkos::parallel_for(
145 "radial profiles avg", num_global_shells, KOKKOS_LAMBDA(
int r ) {
185template <
typename ScalarType >
194 throw std::runtime_error(
"Radial profiles and radii do not have the same number of shells." );
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_ );
209 for (
int r = 0; r < radii.size(); r++ )
212 { {
"tag",
"radial_profiles" },
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 ) } } );
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,
246 const std::vector< ProfileInRadiiType >& profile_radii,
247 const std::vector< ProfileInValueType >& profile_values )
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 );
253 auto coords_radii_host = Kokkos::create_mirror_view_and_copy( Kokkos::HostSpace{}, coords_radii );
257 std::vector< int > idx( profile_radii.size() );
258 std::iota( idx.begin(), idx.end(), 0 );
261 std::sort( idx.begin(), idx.end(), [&](
int a,
int b ) { return profile_radii[a] < profile_radii[b]; } );
264 std::vector< ProfileInRadiiType > profile_radii_sorted( profile_radii );
265 std::vector< ProfileInValueType > profile_values_sorted( profile_values );
267 for (
size_t i = 0; i < idx.size(); i++ )
269 profile_radii_sorted[i] = profile_radii[idx[i]];
270 profile_values_sorted[i] = profile_values[idx[i]];
273 for (
int local_subdomain_id = 0; local_subdomain_id < profile_data_host.extent( 0 ); local_subdomain_id++ )
275 const size_t subdomain_size_r = profile_data_host.extent( 1 );
276 const size_t input_profile_size = profile_radii_sorted.size();
278 if ( input_profile_size == 0 )
283 if ( input_profile_size == 1 )
285 for (
int r = 0; r < subdomain_size_r; r++ )
287 profile_data_host( local_subdomain_id, r ) = profile_values_sorted[0];
294 for (
size_t r = 0; r < subdomain_size_r; ++r )
296 double x = coords_radii_host( local_subdomain_id, r );
299 while ( i + 1 < input_profile_size && profile_radii_sorted[i + 1] < x )
305 if ( x <= profile_radii_sorted.front() )
307 profile_data_host( local_subdomain_id, r ) = profile_values_sorted.front();
312 if ( x >= profile_radii_sorted.back() )
314 profile_data_host( local_subdomain_id, r ) = profile_values_sorted.back();
320 profile_radii_sorted[i],
321 profile_radii_sorted[i + 1],
322 profile_values_sorted[i],
323 profile_values_sorted[i + 1],
329 Kokkos::deep_copy( profile_data, profile_data_host );
348template < std::
floating_po
int Gr
idDataType, std::
floating_po
int ProfileOutDataType =
double >
350 const std::string& filename,
351 const std::string& key_radii,
352 const std::string& key_values,
356 if ( profile_table_result.is_err() )
359 Kokkos::abort(
"Error reading csv file into table." );
361 const auto& profile_table = profile_table_result.unwrap();
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 );
367 "radial_profile_" + key_values, coords_radii, profile_radii, profile_values );
376template < std::
floating_po
int Gr
idDataType, std::
floating_po
int ProfileOutDataType =
double >
379 const ProfileOutDataType& value )
382 "radial_profile_constant", coords_radii.extent( 0 ), coords_radii.extent( 1 ) );
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