|
| template<typename ScalarType > |
| 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!).
|
| |
| template<typename ScalarType > |
| 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.
|
| |
| template<std::floating_point GridDataType, std::floating_point ProfileInRadiiType, std::floating_point ProfileInValueType, std::floating_point ProfileOutDataType = double> |
| 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 (local_subdomain_id, subdomain_size_r) for straightforward use in kernels.
|
| |
| template<std::floating_point GridDataType, std::floating_point ProfileOutDataType = double> |
| 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.
|
| |
| template<std::floating_point GridDataType, std::floating_point ProfileOutDataType = double> |
| grid::Grid2DDataScalar< ProfileOutDataType > | interpolate_constant_radial_profile (const grid::Grid2DDataScalar< GridDataType > &coords_radii, const ProfileOutDataType &value) |
| | Interpolating a constant radial profile into a grid.
|
| |
| template<std::floating_point ScalarTypeCoeff, std::floating_point ScalarTypeGrid> |
| grid::Grid3DDataScalar< ScalarTypeCoeff > | spherical_harmonics_coefficients_grid (const int degree_l, const int order_m, const grid::Grid3DDataVec< ScalarTypeGrid, 3 > &coords_shell) |
| | Interpolates spherical harmonics coefficients of a specific order and degree into a grid (on the host), then copies it to the device, then returns the device view.
|
| |
template<std::floating_point GridDataType, std::floating_point ProfileInRadiiType, std::floating_point ProfileInValueType, std::floating_point ProfileOutDataType = double>
| grid::Grid2DDataScalar< ProfileOutDataType > terra::shell::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 (local_subdomain_id, subdomain_size_r) for straightforward use in kernels.
This function takes care of handling possibly duplicated nodes at subdomain boundaries.
- Note
- This will clamp values outside the passed radial profile values to the first and last values in the vector respectively.
- Parameters
-
| profile_out_label | label of the returned Kokkos::View |
| coords_radii | radii of each node for all local subdomains - see terra::grid::shell::subdomain_shell_radii(), the output grid data will have the same extents |
| profile_radii | input profile: vector of radii |
| profile_values | input profile: vector of values to interpolate |
- Returns
- Kokkos::View with the same dimensions of coords_radii, populated with interpolated values per subdomain
template<typename ScalarType >
Compute radial profiles (min, max, sum, count, average) for a Q1 scalar field (on device!).
Computes the radial profiles for a Q1 scalar field by iterating over all radial shells. The profiles include minimum, maximum, sum, and count of nodes in each shell. The radial shells are defined by the radial dimension of the Q1 scalar field. The output is a RadialProfiles struct with device-side arrays of size num_shells which contain:
- Minimum value in the shell
- Maximum value in the shell
- Sum of values in the shell
- Count of nodes in the shell
- Average of nodes in the shell (sum / count) Performs reduction per process on the device and also inter-device reduction with MPI_Allreduce. All processes will carry the same result.
Mask data is used internally to filter out non-owned nodes. So this really only reduces owned nodes.
- Note
- The returned Views are still on the device. To convert it to a util::Table for output, use the
radial_profiles_to_table() function. So unless you really want to compute any further data on the device using the profile, you like always want to call the radial_profiles_to_table() function as outlined below.
To nicely format the output, use the radial_profiles_to_table() function, e.g., via
auto radii = domain_info.radii();
std::ofstream out( "radial_profiles.csv" );
table.print_csv( out );
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:195
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
- Template Parameters
-
| ScalarType | Scalar type of the field. |
- Parameters
-
| data | Q1 scalar field data. |
| subdomain_shell_idx | global shell indices array with layout subdomain_shell_idx( local_subdomain_id, node_r_idx ) = global_shell_idx
compute with terra::grid::shell::subdomain_shell_idx. |
| num_global_shells | number of global shells |
- Returns
- RadialProfiles struct containing [min, max, sum, count, average] for each radial shell (still on the device).
template<typename ScalarType >
| util::Table terra::shell::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.
Converts the radial profile data (min, max, sum, avg, count) into a util::Table.
This table can then be used for further analysis or output to CSV/JSON. The table will have the following columns:
- tag: "radial_profiles"
- shell_idx: Index of the radial shell
- radius: Radius of the shell
- min: Minimum value in the shell
- max: Maximum value in the shell
- sum: Average value in the shell
- avg: Average value in the shell
- cnt: Count of nodes in the shell
To use this function, you can compute the radial profiles using radial_profiles() and then convert the result to a table using this function:
auto radii = domain_info.radii();
std::ofstream out( "radial_profiles.csv" );
table.print_csv( out );
- Template Parameters
-
| ScalarType | Scalar type of the profile data. |
- Parameters
-
| radial_profiles | RadialProfiles struct containing radial profile statistics. Compute this using radial_profiles() function. Data is expected to be on the device still. It is copied to host for table creation in this function. |
| radii | Vector of shell radii. Can for instance be obtained from the DomainInfo. |
- Returns
- Table with columns: tag, shell_idx, radius, min, max, avg, cnt.