Thermodynamics Observer¶
Header file: <libs/observers/thermo_observer.hpp>
[source]
-
template<typename Store, typename FunctorFunc>
CollectDataForDataset<Store> auto CollectThermoVariable(const Dataset<Store> &dataset, const FunctorFunc ffunc, const std::string_view name, const std::string_view units, const double scale_factor, const size_t maxchunk, const size_t ngbxs)¶ Constructs type sastifying the CollectDataForDataset concept for a given Store (using an instance of the GenericCollectData class) which writes a thermodynamic variable to an Xarray in a dataset.
Function return type writes a thermodyanmic varaible “name” to an Xarray as a 4-byte floating point type by collecting data according to the given FunctorFunc from within a Kokkos::parallel_for loop over gridboxes with a range policy.
- Parameters:
dataset – The dataset to write the variable to.
ffunc – The functor function to collect the variable from within a parallel range policy over gridboxes.
name – The name of the new coordinate.
units – The units of the coordinate.
scale_factor – The scale factor of the coordinate data.
maxchunk – The maximum chunk size (number of elements).
ngbxs – The number of gridboxes.
- Returns:
CollectDataForDataset<Store> An instance satisfying the CollectDataForDataset concept for collecting a 2-D floating point variable (e.g. a thermodynamic variable) from each gridbox.
-
struct PressFunc¶
Functor operator to perform a copy of the pressure from the state of each gridbox to d_data within Kokkos::parallel_for loop over gridboxes with range policy.
Signature of operator such that type can be used by GenericCollectData struct for FunctorFunc.
Note: Conversion of press from double (8 bytes) to single precision float (4 bytes).
- Param ii:
The index of the gridbox.
- Param d_gbxs:
The view of gridboxes on device.
- Param totsupers:
The view of superdroplets on device.
- Param d_data:
The mirror view buffer for the pressure in each gridbox.
Public Functions
-
inline void operator()(const size_t ii, viewd_constgbx d_gbxs, const viewd_constsupers totsupers, Buffer<float>::mirrorviewd_buffer d_data) const¶
-
struct TempFunc¶
Functor operator to perform a copy of the temperature from the state of each gridbox to d_data within Kokkos::parallel_for loop over gridboxes with range policy.
Signature of operator such that type can be used by GenericCollectData struct for FunctorFunc.
Note: Conversion of temp from double (8 bytes) to single precision float (4 bytes).
- Param ii:
The index of the gridbox.
- Param d_gbxs:
The view of gridboxes on device.
- Param totsupers:
The view of superdroplets on device.
- Param d_data:
The mirror view buffer for the temperature in each gridbox.
Public Functions
-
inline void operator()(const size_t ii, viewd_constgbx d_gbxs, const viewd_constsupers totsupers, Buffer<float>::mirrorviewd_buffer d_data) const¶
-
struct QvapFunc¶
Functor operator to perform a copy of the vapour mass mixing ratio “qvap” from the state of each gridbox to d_data within Kokkos::parallel_for loop over gridboxes with range policy.
Signature of operator such that type can be used by GenericCollectData struct for FunctorFunc.
Note: Conversion of qvap from double (8 bytes) to single precision float (4 bytes).
- Param ii:
The index of the gridbox.
- Param d_gbxs:
The view of gridboxes on device.
- Param totsupers:
The view of superdroplets on device.
- Param d_data:
The mirror view buffer for qvap from each gridbox.
Public Functions
-
inline void operator()(const size_t ii, viewd_constgbx d_gbxs, const viewd_constsupers totsupers, Buffer<float>::mirrorviewd_buffer d_data) const¶
-
struct QcondFunc¶
Functor operator to perform a copy of the liquid mass mixing ratio “qcond” from the state of each gridbox to d_data within Kokkos::parallel_for loop over gridboxes with range policy.
Signature of operator such that type can be used by GenericCollectData struct for FunctorFunc.
Note: Conversion of qcond from double (8 bytes) to single precision float (4 bytes).
- Param ii:
The index of the gridbox.
- Param d_gbxs:
The view of gridboxes on device.
- Param totsupers:
The view of superdroplets on device.
- Param d_data:
The mirror view buffer for qcond from each gridbox.
Public Functions
-
inline void operator()(const size_t ii, viewd_constgbx d_gbxs, const viewd_constsupers totsupers, Buffer<float>::mirrorviewd_buffer d_data) const¶
-
template<typename Store>
inline CollectDataForDataset<Store> auto CollectThermo(const Dataset<Store> &dataset, const size_t maxchunk, const size_t ngbxs)¶ Constructs a type satisyfing the CollectDataForDataset concept for collecting multiple thermodyanmcis variables from each gridbox and writing them to a dataset.
This function combines CollectDataForDataset types for many thermodynamic variables from each gridbox (e.g. press, temp, qvap, qcond, etc.) using instances of the GenericCollectData class.
- Template Parameters:
Store – The type of the dataset store.
- Parameters:
dataset – The dataset to write the wind velocity components to.
maxchunk – The maximum chunk size (number of elements).
ngbxs – The number of gridboxes.
- Returns:
CollectDataForDataset<Store> An instance of CollectDataForDataset for collecting thermodynamics from the state of each gridbox.
-
template<typename Store>
inline Observer auto ThermoObserver(const unsigned int interval, const Dataset<Store> &dataset, const size_t maxchunk, const size_t ngbxs)¶ Constructs an observer which writes thermodyanmcis from each gridbox (e.g. press, temp, qvap, etc.) at start of each observation timestep to a arrays with a constant observation timestep “interval”.
- Template Parameters:
Store – Type of store for dataset.
- Parameters:
interval – Observation timestep.
dataset – Dataset to write time data to.
maxchunk – Maximum number of elements in a chunk (1-D vector size).
ngbxs – The number of gridboxes.
- Returns:
Observer An observer instance for writing thermodynamic variables from each gridbox.